20 #ifndef EIGEN_BDCSVD_H 21 #define EIGEN_BDCSVD_H 27 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE 28 IOFormat bdcsvdfmt(8, 0,
", ",
"\n",
" [",
"]");
31 template<
typename _MatrixType>
class BDCSVD;
35 template<
typename _MatrixType>
38 typedef _MatrixType MatrixType;
66 template<
typename _MatrixType>
77 typedef _MatrixType MatrixType;
78 typedef typename MatrixType::Scalar Scalar;
81 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
82 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
83 DiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_DYNAMIC(RowsAtCompileTime, ColsAtCompileTime),
84 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
85 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
86 MaxDiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(MaxRowsAtCompileTime, MaxColsAtCompileTime),
87 MatrixOptions = MatrixType::Options
90 typedef typename Base::MatrixUType MatrixUType;
91 typedef typename Base::MatrixVType MatrixVType;
92 typedef typename Base::SingularValuesType SingularValuesType;
118 : m_algoswap(16), m_numIters(0)
120 allocate(rows, cols, computationOptions);
133 BDCSVD(
const MatrixType& matrix,
unsigned int computationOptions = 0)
134 : m_algoswap(16), m_numIters(0)
136 compute(matrix, computationOptions);
153 BDCSVD& compute(
const MatrixType& matrix,
unsigned int computationOptions);
163 return compute(matrix, this->m_computationOptions);
166 void setSwitchSize(
int s)
168 eigen_assert(s>3 &&
"BDCSVD the size of the algo switch has to be greater than 3");
173 void allocate(
Index rows,
Index cols,
unsigned int computationOptions);
175 void computeSVDofM(
Index firstCol,
Index n, MatrixXr& U, VectorType& singVals, MatrixXr& V);
176 void computeSingVals(
const ArrayRef& col0,
const ArrayRef& diag,
const IndicesRef& perm, VectorType& singVals, ArrayRef shifts, ArrayRef mus);
177 void perturbCol0(
const ArrayRef& col0,
const ArrayRef& diag,
const IndicesRef& perm,
const VectorType& singVals,
const ArrayRef& shifts,
const ArrayRef& mus, ArrayRef zhat);
178 void computeSingVecs(
const ArrayRef& zhat,
const ArrayRef& diag,
const IndicesRef& perm,
const VectorType& singVals,
const ArrayRef& shifts,
const ArrayRef& mus, MatrixXr& U, MatrixXr& V);
182 template<
typename HouseholderU,
typename HouseholderV,
typename NaiveU,
typename NaiveV>
183 void copyUV(
const HouseholderU &householderU,
const HouseholderV &householderV,
const NaiveU &naiveU,
const NaiveV &naivev);
185 static RealScalar secularEq(RealScalar x,
const ArrayRef& col0,
const ArrayRef& diag,
const IndicesRef &perm,
const ArrayRef& diagShifted, RealScalar shift);
188 MatrixXr m_naiveU, m_naiveV;
192 ArrayXi m_workspaceI;
194 bool m_isTranspose, m_compU, m_compV;
196 using Base::m_singularValues;
197 using Base::m_diagSize;
198 using Base::m_computeFullU;
199 using Base::m_computeFullV;
200 using Base::m_computeThinU;
201 using Base::m_computeThinV;
202 using Base::m_matrixU;
203 using Base::m_matrixV;
204 using Base::m_isInitialized;
205 using Base::m_nonzeroSingularValues;
213 template<
typename MatrixType>
216 m_isTranspose = (cols > rows);
218 if (Base::allocate(rows, cols, computationOptions))
221 m_computed = MatrixXr::Zero(m_diagSize + 1, m_diagSize );
222 m_compU = computeV();
223 m_compV = computeU();
225 std::swap(m_compU, m_compV);
227 if (m_compU) m_naiveU = MatrixXr::Zero(m_diagSize + 1, m_diagSize + 1 );
228 else m_naiveU = MatrixXr::Zero(2, m_diagSize + 1 );
230 if (m_compV) m_naiveV = MatrixXr::Zero(m_diagSize, m_diagSize);
232 m_workspace.resize((m_diagSize+1)*(m_diagSize+1)*3);
233 m_workspaceI.resize(3*m_diagSize);
236 template<
typename MatrixType>
239 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE 240 std::cout <<
"\n\n\n======================================================================================================================\n\n\n";
242 allocate(matrix.rows(), matrix.cols(), computationOptions);
245 const RealScalar considerZero = (std::numeric_limits<RealScalar>::min)();
248 if(matrix.cols() < m_algoswap)
252 if(computeU()) m_matrixU = jsvd.
matrixU();
253 if(computeV()) m_matrixV = jsvd.
matrixV();
256 m_isInitialized =
true;
261 RealScalar scale = matrix.cwiseAbs().maxCoeff();
262 if(scale==RealScalar(0)) scale = RealScalar(1);
264 if (m_isTranspose) copy = matrix.adjoint()/scale;
265 else copy = matrix/scale;
275 m_computed.topRows(m_diagSize) = bid.bidiagonal().toDenseMatrix().transpose();
276 m_computed.template bottomRows<1>().setZero();
277 divide(0, m_diagSize - 1, 0, 0, 0);
280 for (
int i=0; i<m_diagSize; i++)
282 RealScalar a = abs(m_computed.coeff(i, i));
283 m_singularValues.coeffRef(i) = a * scale;
286 m_nonzeroSingularValues = i;
287 m_singularValues.tail(m_diagSize - i - 1).setZero();
290 else if (i == m_diagSize - 1)
292 m_nonzeroSingularValues = i + 1;
297 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE 301 if(m_isTranspose) copyUV(bid.householderV(), bid.householderU(), m_naiveV, m_naiveU);
302 else copyUV(bid.householderU(), bid.householderV(), m_naiveU, m_naiveV);
304 m_isInitialized =
true;
309 template<
typename MatrixType>
310 template<
typename HouseholderU,
typename HouseholderV,
typename NaiveU,
typename NaiveV>
311 void BDCSVD<MatrixType>::copyUV(
const HouseholderU &householderU,
const HouseholderV &householderV,
const NaiveU &naiveU,
const NaiveV &naiveV)
316 Index Ucols = m_computeThinU ? m_diagSize : householderU.cols();
317 m_matrixU = MatrixX::Identity(householderU.cols(), Ucols);
318 m_matrixU.topLeftCorner(m_diagSize, m_diagSize) = naiveV.template cast<Scalar>().topLeftCorner(m_diagSize, m_diagSize);
319 householderU.applyThisOnTheLeft(m_matrixU);
323 Index Vcols = m_computeThinV ? m_diagSize : householderV.cols();
324 m_matrixV = MatrixX::Identity(householderV.cols(), Vcols);
325 m_matrixV.topLeftCorner(m_diagSize, m_diagSize) = naiveU.template cast<Scalar>().topLeftCorner(m_diagSize, m_diagSize);
326 householderV.applyThisOnTheLeft(m_matrixV);
338 template<
typename MatrixType>
352 for(
Index j=0; j<n; ++j)
354 if( (A.col(j).head(n1).array()!=0).any() )
356 A1.col(k1) = A.col(j).head(n1);
357 B1.row(k1) = B.row(j);
360 if( (A.col(j).tail(n2).array()!=0).any() )
362 A2.col(k2) = A.col(j).tail(n2);
363 B2.row(k2) = B.row(j);
368 A.topRows(n1).noalias() = A1.leftCols(k1) * B1.topRows(k1);
369 A.bottomRows(n2).noalias() = A2.leftCols(k2) * B2.topRows(k2);
389 template<
typename MatrixType>
396 const Index n = lastCol - firstCol + 1;
398 const RealScalar considerZero = (std::numeric_limits<RealScalar>::min)();
402 RealScalar lambda, phi, c0, s0;
411 m_naiveU.block(firstCol, firstCol, n + 1, n + 1).real() = b.
matrixU();
414 m_naiveU.row(0).segment(firstCol, n + 1).real() = b.
matrixU().row(0);
415 m_naiveU.row(1).segment(firstCol, n + 1).real() = b.
matrixU().row(n);
417 if (m_compV) m_naiveV.block(firstRowW, firstColW, n, n).real() = b.
matrixV();
418 m_computed.block(firstCol + shift, firstCol + shift, n + 1, n).setZero();
419 m_computed.diagonal().segment(firstCol + shift, n) = b.
singularValues().head(n);
423 alphaK = m_computed(firstCol + k, firstCol + k);
424 betaK = m_computed(firstCol + k + 1, firstCol + k);
428 divide(k + 1 + firstCol, lastCol, k + 1 + firstRowW, k + 1 + firstColW, shift);
429 divide(firstCol, k - 1 + firstCol, firstRowW, firstColW + 1, shift + 1);
433 lambda = m_naiveU(firstCol + k, firstCol + k);
434 phi = m_naiveU(firstCol + k + 1, lastCol + 1);
438 lambda = m_naiveU(1, firstCol + k);
439 phi = m_naiveU(0, lastCol + 1);
441 r0 = sqrt((abs(alphaK * lambda) * abs(alphaK * lambda)) + abs(betaK * phi) * abs(betaK * phi));
444 l = m_naiveU.row(firstCol + k).segment(firstCol, k);
445 f = m_naiveU.row(firstCol + k + 1).segment(firstCol + k + 1, n - k - 1);
449 l = m_naiveU.row(1).segment(firstCol, k);
450 f = m_naiveU.row(0).segment(firstCol + k + 1, n - k - 1);
452 if (m_compV) m_naiveV(firstRowW+k, firstColW) = 1;
460 c0 = alphaK * lambda / r0;
461 s0 = betaK * phi / r0;
464 #ifdef EIGEN_BDCSVD_SANITY_CHECKS 465 assert(m_naiveU.allFinite());
466 assert(m_naiveV.allFinite());
467 assert(m_computed.allFinite());
472 MatrixXr q1 (m_naiveU.col(firstCol + k).segment(firstCol, k + 1));
474 for (
Index i = firstCol + k - 1; i >= firstCol; i--)
475 m_naiveU.col(i + 1).segment(firstCol, k + 1) = m_naiveU.col(i).segment(firstCol, k + 1);
477 m_naiveU.col(firstCol).segment( firstCol, k + 1) = (q1 * c0);
479 m_naiveU.col(lastCol + 1).segment(firstCol, k + 1) = (q1 * ( - s0));
481 m_naiveU.col(firstCol).segment(firstCol + k + 1, n - k) = m_naiveU.col(lastCol + 1).segment(firstCol + k + 1, n - k) * s0;
483 m_naiveU.col(lastCol + 1).segment(firstCol + k + 1, n - k) *= c0;
487 RealScalar q1 = m_naiveU(0, firstCol + k);
489 for (
Index i = firstCol + k - 1; i >= firstCol; i--)
490 m_naiveU(0, i + 1) = m_naiveU(0, i);
492 m_naiveU(0, firstCol) = (q1 * c0);
494 m_naiveU(0, lastCol + 1) = (q1 * ( - s0));
496 m_naiveU(1, firstCol) = m_naiveU(1, lastCol + 1) *s0;
498 m_naiveU(1, lastCol + 1) *= c0;
499 m_naiveU.row(1).segment(firstCol + 1, k).setZero();
500 m_naiveU.row(0).segment(firstCol + k + 1, n - k - 1).setZero();
503 #ifdef EIGEN_BDCSVD_SANITY_CHECKS 504 assert(m_naiveU.allFinite());
505 assert(m_naiveV.allFinite());
506 assert(m_computed.allFinite());
509 m_computed(firstCol + shift, firstCol + shift) = r0;
510 m_computed.col(firstCol + shift).segment(firstCol + shift + 1, k) = alphaK * l.transpose().real();
511 m_computed.col(firstCol + shift).segment(firstCol + shift + k + 1, n - k - 1) = betaK * f.transpose().real();
513 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE 514 ArrayXr tmp1 = (m_computed.block(firstCol+shift, firstCol+shift, n, n)).jacobiSvd().singularValues();
517 deflation(firstCol, lastCol, k, firstRowW, firstColW, shift);
518 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE 519 ArrayXr tmp2 = (m_computed.block(firstCol+shift, firstCol+shift, n, n)).jacobiSvd().singularValues();
520 std::cout <<
"\n\nj1 = " << tmp1.transpose().format(bdcsvdfmt) <<
"\n";
521 std::cout <<
"j2 = " << tmp2.transpose().format(bdcsvdfmt) <<
"\n\n";
522 std::cout <<
"err: " << ((tmp1-tmp2).abs()>1e-12*tmp2.abs()).transpose() <<
"\n";
523 static int count = 0;
524 std::cout <<
"# " << ++count <<
"\n\n";
525 assert((tmp1-tmp2).matrix().norm() < 1e-14*tmp2.matrix().norm());
533 computeSVDofM(firstCol + shift, n, UofSVD, singVals, VofSVD);
535 #ifdef EIGEN_BDCSVD_SANITY_CHECKS 536 assert(UofSVD.allFinite());
537 assert(VofSVD.allFinite());
541 structured_update(m_naiveU.block(firstCol, firstCol, n + 1, n + 1), UofSVD, (n+2)/2);
545 tmp.noalias() = m_naiveU.middleCols(firstCol, n+1) * UofSVD;
546 m_naiveU.middleCols(firstCol, n + 1) = tmp;
549 if (m_compV) structured_update(m_naiveV.block(firstRowW, firstColW, n, n), VofSVD, (n+1)/2);
551 #ifdef EIGEN_BDCSVD_SANITY_CHECKS 552 assert(m_naiveU.allFinite());
553 assert(m_naiveV.allFinite());
554 assert(m_computed.allFinite());
557 m_computed.block(firstCol + shift, firstCol + shift, n, n).setZero();
558 m_computed.block(firstCol + shift, firstCol + shift, n, n).diagonal() = singVals;
569 template <
typename MatrixType>
572 const RealScalar considerZero = (std::numeric_limits<RealScalar>::min)();
574 ArrayRef col0 = m_computed.col(firstCol).segment(firstCol, n);
575 m_workspace.head(n) = m_computed.block(firstCol, firstCol, n, n).diagonal();
576 ArrayRef diag = m_workspace.head(n);
582 if (m_compV) V.
resize(n, n);
584 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE 585 if (col0.hasNaN() || diag.hasNaN())
586 std::cout <<
"\n\nHAS NAN\n\n";
593 while(actual_n>1 && diag(actual_n-1)==0) --actual_n;
595 for(
Index k=0;k<actual_n;++k)
596 if(abs(col0(k))>considerZero)
597 m_workspaceI(m++) = k;
604 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE 605 std::cout <<
"computeSVDofM using:\n";
606 std::cout <<
" z: " << col0.transpose() <<
"\n";
607 std::cout <<
" d: " << diag.transpose() <<
"\n";
611 computeSingVals(col0, diag, perm, singVals, shifts, mus);
613 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE 614 std::cout <<
" j: " << (m_computed.block(firstCol, firstCol, n, n)).jacobiSvd().singularValues().transpose().reverse() <<
"\n\n";
615 std::cout <<
" sing-val: " << singVals.transpose() <<
"\n";
616 std::cout <<
" mu: " << mus.transpose() <<
"\n";
617 std::cout <<
" shift: " << shifts.transpose() <<
"\n";
621 while(actual_n>1 && abs(col0(actual_n-1))<considerZero) --actual_n;
622 std::cout <<
"\n\n mus: " << mus.head(actual_n).transpose() <<
"\n\n";
623 std::cout <<
" check1 (expect0) : " << ((singVals.array()-(shifts+mus)) / singVals.array()).head(actual_n).transpose() <<
"\n\n";
624 std::cout <<
" check2 (>0) : " << ((singVals.array()-diag) / singVals.array()).head(actual_n).transpose() <<
"\n\n";
625 std::cout <<
" check3 (>0) : " << ((diag.segment(1,actual_n-1)-singVals.head(actual_n-1).array()) / singVals.head(actual_n-1).array()).transpose() <<
"\n\n\n";
626 std::cout <<
" check4 (>0) : " << ((singVals.segment(1,actual_n-1)-singVals.head(actual_n-1))).transpose() <<
"\n\n\n";
630 #ifdef EIGEN_BDCSVD_SANITY_CHECKS 631 assert(singVals.allFinite());
632 assert(mus.allFinite());
633 assert(shifts.allFinite());
637 perturbCol0(col0, diag, perm, singVals, shifts, mus, zhat);
638 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE 639 std::cout <<
" zhat: " << zhat.transpose() <<
"\n";
642 #ifdef EIGEN_BDCSVD_SANITY_CHECKS 643 assert(zhat.allFinite());
646 computeSingVecs(zhat, diag, perm, singVals, shifts, mus, U, V);
648 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE 649 std::cout <<
"U^T U: " << (U.transpose() * U -
MatrixXr(MatrixXr::Identity(U.cols(),U.cols()))).norm() <<
"\n";
650 std::cout <<
"V^T V: " << (V.transpose() * V -
MatrixXr(MatrixXr::Identity(V.cols(),V.cols()))).norm() <<
"\n";
653 #ifdef EIGEN_BDCSVD_SANITY_CHECKS 654 assert(U.allFinite());
655 assert(V.allFinite());
656 assert((U.transpose() * U -
MatrixXr(MatrixXr::Identity(U.cols(),U.cols()))).norm() < 1e-14 * n);
657 assert((V.transpose() * V -
MatrixXr(MatrixXr::Identity(V.cols(),V.cols()))).norm() < 1e-14 * n);
658 assert(m_naiveU.allFinite());
659 assert(m_naiveV.allFinite());
660 assert(m_computed.allFinite());
665 for(
Index i=0; i<actual_n-1; ++i)
667 if(singVals(i)>singVals(i+1))
670 swap(singVals(i),singVals(i+1));
671 U.col(i).swap(U.col(i+1));
672 if(m_compV) V.col(i).swap(V.col(i+1));
678 singVals.head(actual_n).reverseInPlace();
679 U.leftCols(actual_n).rowwise().reverseInPlace();
680 if (m_compV) V.leftCols(actual_n).rowwise().reverseInPlace();
682 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE 684 std::cout <<
" * j: " << jsvd.singularValues().transpose() <<
"\n\n";
685 std::cout <<
" * sing-val: " << singVals.transpose() <<
"\n";
690 template <
typename MatrixType>
693 Index m = perm.size();
695 for(
Index i=0; i<m; ++i)
698 res += numext::abs2(col0(j)) / ((diagShifted(j) - mu) * (diag(j) + shift + mu));
704 template <
typename MatrixType>
711 Index n = col0.size();
713 while(actual_n>1 && col0(actual_n-1)==0) --actual_n;
715 for (
Index k = 0; k < n; ++k)
717 if (col0(k) == 0 || actual_n==1)
721 singVals(k) = k==0 ? col0(0) : diag(k);
723 shifts(k) = k==0 ? col0(0) : diag(k);
728 RealScalar left = diag(k);
731 right = (diag(actual_n-1) + col0.matrix().norm());
736 while(col0(l)==0) { ++l; eigen_internal_assert(l<actual_n); }
741 RealScalar mid = left + (right-left) / 2;
742 RealScalar fMid = secularEq(mid, col0, diag, perm, diag, 0);
743 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE 744 std::cout << right-left <<
"\n";
745 std::cout <<
"fMid = " << fMid <<
" " << secularEq(mid-left, col0, diag, perm, diag-left, left) <<
" " << secularEq(mid-right, col0, diag, perm, diag-right, right) <<
"\n";
746 std::cout <<
" = " << secularEq(0.1*(left+right), col0, diag, perm, diag, 0)
747 <<
" " << secularEq(0.2*(left+right), col0, diag, perm, diag, 0)
748 <<
" " << secularEq(0.3*(left+right), col0, diag, perm, diag, 0)
749 <<
" " << secularEq(0.4*(left+right), col0, diag, perm, diag, 0)
750 <<
" " << secularEq(0.49*(left+right), col0, diag, perm, diag, 0)
751 <<
" " << secularEq(0.5*(left+right), col0, diag, perm, diag, 0)
752 <<
" " << secularEq(0.51*(left+right), col0, diag, perm, diag, 0)
753 <<
" " << secularEq(0.6*(left+right), col0, diag, perm, diag, 0)
754 <<
" " << secularEq(0.7*(left+right), col0, diag, perm, diag, 0)
755 <<
" " << secularEq(0.8*(left+right), col0, diag, perm, diag, 0)
756 <<
" " << secularEq(0.9*(left+right), col0, diag, perm, diag, 0) <<
"\n";
758 RealScalar shift = (k == actual_n-1 || fMid > 0) ? left : right;
762 diagShifted = diag - shift;
765 RealScalar muPrev, muCur;
768 muPrev = (right - left) * RealScalar(0.1);
769 if (k == actual_n-1) muCur = right - left;
770 else muCur = (right - left) * RealScalar(0.5);
774 muPrev = -(right - left) * RealScalar(0.1);
775 muCur = -(right - left) * RealScalar(0.5);
778 RealScalar fPrev = secularEq(muPrev, col0, diag, perm, diagShifted, shift);
779 RealScalar fCur = secularEq(muCur, col0, diag, perm, diagShifted, shift);
780 if (abs(fPrev) < abs(fCur))
788 bool useBisection = fPrev*fCur>0;
794 RealScalar a = (fCur - fPrev) / (1/muCur - 1/muPrev);
795 RealScalar b = fCur - a / muCur;
797 RealScalar muZero = -a/b;
798 RealScalar fZero = secularEq(muZero, col0, diag, perm, diagShifted, shift);
806 if (shift == left && (muCur < 0 || muCur > right - left)) useBisection =
true;
807 if (shift == right && (muCur < -(right - left) || muCur > 0)) useBisection =
true;
808 if (abs(fCur)>abs(fPrev)) useBisection =
true;
814 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE 815 std::cout <<
"useBisection for k = " << k <<
", actual_n = " << actual_n <<
"\n";
817 RealScalar leftShifted, rightShifted;
820 leftShifted = (std::numeric_limits<RealScalar>::min)();
823 rightShifted = (k==actual_n-1) ? right : ((right - left) * RealScalar(0.6));
827 leftShifted = -(right - left) * RealScalar(0.6);
828 rightShifted = -(std::numeric_limits<RealScalar>::min)();
831 RealScalar fLeft = secularEq(leftShifted, col0, diag, perm, diagShifted, shift);
833 #if defined EIGEN_INTERNAL_DEBUGGING || defined EIGEN_BDCSVD_DEBUG_VERBOSE 834 RealScalar fRight = secularEq(rightShifted, col0, diag, perm, diagShifted, shift);
837 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE 838 if(!(fLeft * fRight<0))
840 std::cout <<
"fLeft: " << leftShifted <<
" - " << diagShifted.head(10).transpose() <<
"\n ; " << bool(left==shift) <<
" " << (left-shift) <<
"\n";
841 std::cout << k <<
" : " << fLeft <<
" * " << fRight <<
" == " << fLeft * fRight <<
" ; " << left <<
" - " << right <<
" -> " << leftShifted <<
" " << rightShifted <<
" shift=" << shift <<
"\n";
844 eigen_internal_assert(fLeft * fRight < 0);
848 RealScalar midShifted = (leftShifted + rightShifted) / 2;
849 fMid = secularEq(midShifted, col0, diag, perm, diagShifted, shift);
850 if (fLeft * fMid < 0)
852 rightShifted = midShifted;
856 leftShifted = midShifted;
861 muCur = (leftShifted + rightShifted) / 2;
864 singVals[k] = shift + muCur;
878 template <
typename MatrixType>
884 Index n = col0.size();
885 Index m = perm.size();
891 Index last = perm(m-1);
893 for (
Index k = 0; k < n; ++k)
900 RealScalar dk = diag(k);
901 RealScalar prod = (singVals(last) + dk) * (mus(last) + (shifts(last) - dk));
903 for(
Index l = 0; l<m; ++l)
908 Index j = i<k ? i : perm(l-1);
909 prod *= ((singVals(j)+dk) / ((diag(i)+dk))) * ((mus(j)+(shifts(j)-dk)) / ((diag(i)-dk)));
910 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE 911 if(i!=k && std::abs(((singVals(j)+dk)*(mus(j)+(shifts(j)-dk)))/((diag(i)+dk)*(diag(i)-dk)) - 1) > 0.9 )
912 std::cout <<
" " << ((singVals(j)+dk)*(mus(j)+(shifts(j)-dk)))/((diag(i)+dk)*(diag(i)-dk)) <<
" == (" << (singVals(j)+dk) <<
" * " << (mus(j)+(shifts(j)-dk))
913 <<
") / (" << (diag(i)+dk) <<
" * " << (diag(i)-dk) <<
")\n";
917 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE 918 std::cout <<
"zhat(" << k <<
") = sqrt( " << prod <<
") ; " << (singVals(last) + dk) <<
" * " << mus(last) + shifts(last) <<
" - " << dk <<
"\n";
920 RealScalar tmp = sqrt(prod);
921 zhat(k) = col0(k) > 0 ? tmp : -tmp;
927 template <
typename MatrixType>
932 Index n = zhat.size();
933 Index m = perm.size();
935 for (
Index k = 0; k < n; ++k)
939 U.col(k) = VectorType::Unit(n+1, k);
940 if (m_compV) V.col(k) = VectorType::Unit(n, k);
945 for(
Index l=0;l<m;++l)
948 U(i,k) = zhat(i)/(((diag(i) - shifts(k)) - mus(k)) )/( (diag(i) + singVals[k]));
951 U.col(k).normalize();
956 for(
Index l=1;l<m;++l)
959 V(i,k) = diag(i) * zhat(i) / (((diag(i) - shifts(k)) - mus(k)) )/( (diag(i) + singVals[k]));
962 V.col(k).normalize();
966 U.col(n) = VectorType::Unit(n+1, n);
973 template <
typename MatrixType>
979 Index start = firstCol + shift;
980 RealScalar c = m_computed(start, start);
981 RealScalar s = m_computed(start+i, start);
982 RealScalar r = sqrt(numext::abs2(c) + numext::abs2(s));
985 m_computed(start+i, start+i) = 0;
988 m_computed(start,start) = r;
989 m_computed(start+i, start) = 0;
990 m_computed(start+i, start+i) = 0;
993 if (m_compU) m_naiveU.middleRows(firstCol, size+1).applyOnTheRight(firstCol, firstCol+i, J);
994 else m_naiveU.applyOnTheRight(firstCol, firstCol+i, J);
1002 template <
typename MatrixType>
1009 RealScalar c = m_computed(firstColm+i, firstColm);
1010 RealScalar s = m_computed(firstColm+j, firstColm);
1011 RealScalar r = sqrt(numext::abs2(c) + numext::abs2(s));
1012 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE 1013 std::cout <<
"deflation 4.4: " << i <<
"," << j <<
" -> " << c <<
" " << s <<
" " << r <<
" ; " 1014 << m_computed(firstColm + i-1, firstColm) <<
" " 1015 << m_computed(firstColm + i, firstColm) <<
" " 1016 << m_computed(firstColm + i+1, firstColm) <<
" " 1017 << m_computed(firstColm + i+2, firstColm) <<
"\n";
1018 std::cout << m_computed(firstColm + i-1, firstColm + i-1) <<
" " 1019 << m_computed(firstColm + i, firstColm+i) <<
" " 1020 << m_computed(firstColm + i+1, firstColm+i+1) <<
" " 1021 << m_computed(firstColm + i+2, firstColm+i+2) <<
"\n";
1025 m_computed(firstColm + i, firstColm + i) = m_computed(firstColm + j, firstColm + j);
1030 m_computed(firstColm + i, firstColm) = r;
1031 m_computed(firstColm + j, firstColm + j) = m_computed(firstColm + i, firstColm + i);
1032 m_computed(firstColm + j, firstColm) = 0;
1035 if (m_compU) m_naiveU.middleRows(firstColu, size+1).applyOnTheRight(firstColu + i, firstColu + j, J);
1036 else m_naiveU.applyOnTheRight(firstColu+i, firstColu+j, J);
1037 if (m_compV) m_naiveV.middleRows(firstRowW, size).applyOnTheRight(firstColW + i, firstColW + j, J);
1042 template <
typename MatrixType>
1047 const Index length = lastCol + 1 - firstCol;
1053 const RealScalar considerZero = (std::numeric_limits<RealScalar>::min)();
1054 RealScalar maxDiag = diag.tail((std::max)(
Index(1),length-1)).cwiseAbs().maxCoeff();
1058 #ifdef EIGEN_BDCSVD_SANITY_CHECKS 1059 assert(m_naiveU.allFinite());
1060 assert(m_naiveV.allFinite());
1061 assert(m_computed.allFinite());
1064 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE 1065 std::cout <<
"\ndeflate:" << diag.head(k+1).transpose() <<
" | " << diag.segment(k+1,length-k-1).transpose() <<
"\n";
1069 if (diag(0) < epsilon_coarse)
1071 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE 1072 std::cout <<
"deflation 4.1, because " << diag(0) <<
" < " << epsilon_coarse <<
"\n";
1074 diag(0) = epsilon_coarse;
1078 for (
Index i=1;i<length;++i)
1079 if (abs(col0(i)) < epsilon_strict)
1081 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE 1082 std::cout <<
"deflation 4.2, set z(" << i <<
") to zero because " << abs(col0(i)) <<
" < " << epsilon_strict <<
" (diag(" << i <<
")=" << diag(i) <<
")\n";
1088 for (
Index i=1;i<length; i++)
1089 if (diag(i) < epsilon_coarse)
1091 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE 1092 std::cout <<
"deflation 4.3, cancel z(" << i <<
")=" << col0(i) <<
" because diag(" << i <<
")=" << diag(i) <<
" < " << epsilon_coarse <<
"\n";
1094 deflation43(firstCol, shift, i, length);
1097 #ifdef EIGEN_BDCSVD_SANITY_CHECKS 1098 assert(m_naiveU.allFinite());
1099 assert(m_naiveV.allFinite());
1100 assert(m_computed.allFinite());
1102 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE 1103 std::cout <<
"to be sorted: " << diag.transpose() <<
"\n\n";
1108 bool total_deflation = (col0.tail(length-1).array()<considerZero).all();
1112 Index *permutation = m_workspaceI.data();
1118 for(
Index i=1; i<length; ++i)
1119 if(abs(diag(i))<considerZero)
1120 permutation[p++] = i;
1123 for( ; p < length; ++p)
1125 if (i > k) permutation[p] = j++;
1126 else if (j >= length) permutation[p] = i++;
1127 else if (diag(i) < diag(j)) permutation[p] = j++;
1128 else permutation[p] = i++;
1135 for(
Index i=1; i<length; ++i)
1137 Index pi = permutation[i];
1138 if(abs(diag(pi))<considerZero || diag(0)<diag(pi))
1139 permutation[i-1] = permutation[i];
1142 permutation[i-1] = 0;
1149 Index *realInd = m_workspaceI.data()+length;
1150 Index *realCol = m_workspaceI.data()+2*length;
1152 for(
int pos = 0; pos< length; pos++)
1158 for(
Index i = total_deflation?0:1; i < length; i++)
1160 const Index pi = permutation[length - (total_deflation ? i+1 : i)];
1161 const Index J = realCol[pi];
1165 swap(diag(i), diag(J));
1166 if(i!=0 && J!=0) swap(col0(i), col0(J));
1169 if (m_compU) m_naiveU.col(firstCol+i).segment(firstCol, length + 1).swap(m_naiveU.col(firstCol+J).segment(firstCol, length + 1));
1170 else m_naiveU.col(firstCol+i).segment(0, 2) .swap(m_naiveU.col(firstCol+J).segment(0, 2));
1171 if (m_compV) m_naiveV.col(firstColW + i).segment(firstRowW, length).swap(m_naiveV.col(firstColW + J).segment(firstRowW, length));
1174 const Index realI = realInd[i];
1181 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE 1182 std::cout <<
"sorted: " << diag.transpose().format(bdcsvdfmt) <<
"\n";
1183 std::cout <<
" : " << col0.transpose() <<
"\n\n";
1189 while(i>0 && (abs(diag(i))<considerZero || abs(col0(i))<considerZero)) --i;
1193 #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE 1194 std::cout <<
"deflation 4.4 with i = " << i <<
" because " << (diag(i) - diag(i-1)) <<
" < " <<
NumTraits<RealScalar>::epsilon()*diag(i) <<
"\n";
1196 eigen_internal_assert(abs(diag(i) - diag(i-1))<epsilon_coarse &&
" diagonal entries are not properly sorted");
1197 deflation44(firstCol, firstCol + shift, firstRowW, firstColW, i-1, i, length);
1201 #ifdef EIGEN_BDCSVD_SANITY_CHECKS 1202 for(
Index j=2;j<length;++j)
1203 assert(diag(j-1)<=diag(j) || abs(diag(j))<considerZero);
1206 #ifdef EIGEN_BDCSVD_SANITY_CHECKS 1207 assert(m_naiveU.allFinite());
1208 assert(m_naiveV.allFinite());
1209 assert(m_computed.allFinite());
1220 template<
typename Derived>
BDCSVD & compute(const MatrixType &matrix)
Method performing the decomposition of given matrix using current options.
Definition: BDCSVD.h:161
Used in JacobiSVD to indicate that the square matrix U is to be computed.
Definition: Constants.h:383
BDCSVD< PlainObject > bdcSvd(unsigned int computationOptions=0) const
Definition: BDCSVD.h:1222
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:88
const MatrixUType & matrixU() const
Definition: SVDBase.h:83
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:85
Definition: ForwardDeclarations.h:263
Holds information about the various numeric (i.e.
Definition: NumTraits.h:150
Expression of a fixed-size or dynamic-size sub-vector.
Definition: ForwardDeclarations.h:87
BDCSVD(Index rows, Index cols, unsigned int computationOptions=0)
Default Constructor with memory preallocation.
Definition: BDCSVD.h:117
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
Resizes *this to a rows x cols matrix.
Definition: PlainObjectBase.h:273
Base class of SVD algorithms.
Definition: SVDBase.h:48
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
Definition: Constants.h:235
BDCSVD & compute(const MatrixType &matrix, unsigned int computationOptions)
Method performing the decomposition of given matrix using custom options.
Definition: BDCSVD.h:237
BDCSVD(const MatrixType &matrix, unsigned int computationOptions=0)
Constructor performing the decomposition of given matrix.
Definition: BDCSVD.h:133
Definition: UpperBidiagonalization.h:20
EIGEN_DEVICE_FUNC Derived & setZero(Index size)
Resizes to the given size, and sets all coefficients in this expression to zero.
Definition: CwiseNullaryOp.h:515
class Bidiagonal Divide and Conquer SVD
Definition: ForwardDeclarations.h:259
A matrix or vector expression mapping an existing expression.
Definition: Ref.h:190
Definition: BandTriangularSolver.h:13
Expression of a fixed-size or dynamic-size block.
Definition: Block.h:103
Index nonzeroSingularValues() const
Definition: SVDBase.h:118
BDCSVD()
Default Constructor.
Definition: BDCSVD.h:107
Two-sided Jacobi SVD decomposition of a rectangular matrix.
Definition: ForwardDeclarations.h:258
const SingularValuesType & singularValues() const
Definition: SVDBase.h:111
const MatrixVType & matrixV() const
Definition: SVDBase.h:99
Expression of a diagonal/subdiagonal/superdiagonal in a matrix.
Definition: Diagonal.h:63
const int Dynamic
This value means that a positive quantity (e.g., a size) is not known at compile-time, and that instead the value is stored in some runtime variable.
Definition: Constants.h:21
Used in JacobiSVD to indicate that the square matrix V is to be computed.
Definition: Constants.h:387
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:178
Definition: ForwardDeclarations.h:17