25#ifdef EIGEN_BDCSVD_SANITY_CHECKS
26#undef eigen_internal_assert
27#define eigen_internal_assert(X) assert(X);
31#include "./InternalHeaderCheck.h"
33#ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
39#ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
40IOFormat bdcsvdfmt(8, 0,
", ",
"\n",
" [",
"]");
43template <
typename MatrixType_,
int Options>
48template <
typename MatrixType_,
int Options>
49struct traits<BDCSVD<MatrixType_, Options> > : svd_traits<MatrixType_, Options> {
50 typedef MatrixType_ MatrixType;
84template <
typename MatrixType_,
int Options_>
95 typedef MatrixType_ MatrixType;
96 typedef typename Base::Scalar Scalar;
97 typedef typename Base::RealScalar RealScalar;
98 typedef typename NumTraits<RealScalar>::Literal Literal;
102 QRDecomposition = Options & internal::QRPreconditionerBits,
103 ComputationOptions = Options & internal::ComputationOptionsBits,
104 RowsAtCompileTime = Base::RowsAtCompileTime,
105 ColsAtCompileTime = Base::ColsAtCompileTime,
106 DiagSizeAtCompileTime = Base::DiagSizeAtCompileTime,
107 MaxRowsAtCompileTime = Base::MaxRowsAtCompileTime,
108 MaxColsAtCompileTime = Base::MaxColsAtCompileTime,
109 MaxDiagSizeAtCompileTime = Base::MaxDiagSizeAtCompileTime,
110 MatrixOptions = Base::MatrixOptions
113 typedef typename Base::MatrixUType MatrixUType;
114 typedef typename Base::MatrixVType MatrixVType;
115 typedef typename Base::SingularValuesType SingularValuesType;
130 BDCSVD() : m_algoswap(16), m_isTranspose(false), m_compU(false), m_compV(false), m_numIters(0) {}
138 BDCSVD(Index rows, Index cols) : m_algoswap(16), m_numIters(0) {
139 allocate(rows, cols, internal::get_computation_options(Options));
159 BDCSVD(Index rows, Index cols,
unsigned int computationOptions) : m_algoswap(16), m_numIters(0) {
160 internal::check_svd_options_assertions<MatrixType, Options>(computationOptions, rows, cols);
161 allocate(rows, cols, computationOptions);
169 template <
typename Derived>
171 compute_impl(matrix, internal::get_computation_options(Options));
186 template <
typename Derived>
188 BDCSVD(const
MatrixBase<Derived>& matrix,
unsigned int computationOptions) : m_algoswap(16), m_numIters(0) {
189 internal::check_svd_options_assertions<MatrixType, Options>(computationOptions, matrix.rows(), matrix.cols());
190 compute_impl(matrix, computationOptions);
200 template <
typename Derived>
202 return compute_impl(matrix, m_computationOptions);
214 template <
typename Derived>
217 internal::check_svd_options_assertions<MatrixType, Options>(computationOptions, matrix.rows(), matrix.cols());
218 return compute_impl(matrix, computationOptions);
221 void setSwitchSize(
int s) {
222 eigen_assert(s >= 3 &&
"BDCSVD the size of the algo switch has to be at least 3.");
227 template <
typename Derived>
228 BDCSVD& compute_impl(
const MatrixBase<Derived>& matrix,
unsigned int computationOptions);
229 void divide(Index firstCol, Index lastCol, Index firstRowW, Index firstColW, Index shift);
230 void computeSVDofM(Index firstCol, Index n, MatrixXr& U, VectorType& singVals, MatrixXr& V);
231 void computeSingVals(
const ArrayRef& col0,
const ArrayRef& diag,
const IndicesRef& perm, VectorType& singVals,
232 ArrayRef shifts, ArrayRef mus);
233 void perturbCol0(
const ArrayRef& col0,
const ArrayRef& diag,
const IndicesRef& perm,
const VectorType& singVals,
234 const ArrayRef& shifts,
const ArrayRef& mus, ArrayRef zhat);
235 void computeSingVecs(
const ArrayRef& zhat,
const ArrayRef& diag,
const IndicesRef& perm,
const VectorType& singVals,
236 const ArrayRef& shifts,
const ArrayRef& mus, MatrixXr& U, MatrixXr& V);
237 void deflation43(Index firstCol, Index shift, Index i, Index
size);
238 void deflation44(Index firstColu, Index firstColm, Index firstRowW, Index firstColW, Index i, Index j, Index
size);
239 void deflation(Index firstCol, Index lastCol, Index k, Index firstRowW, Index firstColW, Index shift);
240 template <
typename HouseholderU,
typename HouseholderV,
typename NaiveU,
typename NaiveV>
241 void copyUV(
const HouseholderU& householderU,
const HouseholderV& householderV,
const NaiveU& naiveU,
242 const NaiveV& naivev);
243 void structured_update(Block<MatrixXr, Dynamic, Dynamic> A,
const MatrixXr& B, Index n1);
244 static RealScalar secularEq(RealScalar x,
const ArrayRef& col0,
const ArrayRef& diag,
const IndicesRef& perm,
245 const ArrayRef& diagShifted, RealScalar shift);
246 template <
typename SVDType>
247 void computeBaseCase(SVDType& svd, Index n, Index firstCol, Index firstRowW, Index firstColW, Index shift);
250 void allocate(Index rows, Index cols,
unsigned int computationOptions);
251 MatrixXr m_naiveU, m_naiveV;
255 ArrayXi m_workspaceI;
257 bool m_isTranspose, m_compU, m_compV, m_useQrDecomp;
258 JacobiSVD<MatrixType, ComputationOptions> smallSvd;
259 HouseholderQR<MatrixX> qrDecomp;
260 internal::UpperBidiagonalization<MatrixX> bid;
261 MatrixX copyWorkspace;
262 MatrixX reducedTriangle;
264 using Base::m_computationOptions;
265 using Base::m_computeThinU;
266 using Base::m_computeThinV;
268 using Base::m_isInitialized;
269 using Base::m_matrixU;
270 using Base::m_matrixV;
271 using Base::m_nonzeroSingularValues;
272 using Base::m_singularValues;
279template <
typename MatrixType,
int Options>
280void BDCSVD<MatrixType, Options>::allocate(
Index rows,
Index cols,
unsigned int computationOptions) {
281 if (Base::allocate(rows, cols, computationOptions))
return;
283 if (cols < m_algoswap)
284 smallSvd.allocate(rows, cols, Options == 0 ? computationOptions : internal::get_computation_options(Options));
286 m_computed = MatrixXr::Zero(diagSize() + 1, diagSize());
287 m_compU = computeV();
288 m_compV = computeU();
289 m_isTranspose = (cols > rows);
290 if (m_isTranspose) std::swap(m_compU, m_compV);
296 constexpr Index kMinAspectRatio = 4;
297 constexpr bool disableQrDecomp =
static_cast<int>(QRDecomposition) ==
static_cast<int>(
DisableQRDecomposition);
298 m_useQrDecomp = !disableQrDecomp && ((rows / kMinAspectRatio > cols) || (cols / kMinAspectRatio > rows));
301 reducedTriangle =
MatrixX(diagSize(), diagSize());
304 copyWorkspace =
MatrixX(m_isTranspose ? cols : rows, m_isTranspose ? rows : cols);
305 bid = internal::UpperBidiagonalization<MatrixX>(m_useQrDecomp ? diagSize() : copyWorkspace.rows(),
306 m_useQrDecomp ? diagSize() : copyWorkspace.cols());
309 m_naiveU = MatrixXr::Zero(diagSize() + 1, diagSize() + 1);
311 m_naiveU = MatrixXr::Zero(2, diagSize() + 1);
313 if (m_compV) m_naiveV = MatrixXr::Zero(diagSize(), diagSize());
315 m_workspace.resize((diagSize() + 1) * (diagSize() + 1) * 3);
316 m_workspaceI.resize(3 * diagSize());
319template <
typename MatrixType,
int Options>
320template <
typename Derived>
322 unsigned int computationOptions) {
323 EIGEN_STATIC_ASSERT_SAME_MATRIX_SIZE(Derived,
MatrixType);
324 EIGEN_STATIC_ASSERT((std::is_same<typename Derived::Scalar, typename MatrixType::Scalar>::value),
325 Input matrix must have the same
Scalar type as the
BDCSVD object.);
327#ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
328 std::cout <<
"\n\n\n================================================================================================="
329 "=====================\n\n\n";
333 allocate(matrix.rows(), matrix.cols(), computationOptions);
335 const RealScalar considerZero = (std::numeric_limits<RealScalar>::min)();
338 if (matrix.cols() < m_algoswap) {
339 smallSvd.compute(matrix);
340 m_isInitialized =
true;
341 m_info = smallSvd.info();
343 if (computeU()) m_matrixU = smallSvd.matrixU();
344 if (computeV()) m_matrixV = smallSvd.matrixV();
345 m_singularValues = smallSvd.singularValues();
346 m_nonzeroSingularValues = smallSvd.nonzeroSingularValues();
352 RealScalar scale = matrix.cwiseAbs().template maxCoeff<PropagateNaN>();
353 if (!(numext::isfinite)(scale)) {
354 m_isInitialized =
true;
359 if (numext::is_exactly_zero(scale)) scale =
Literal(1);
362 copyWorkspace = matrix.adjoint() / scale;
364 copyWorkspace = matrix / scale;
371 qrDecomp.compute(copyWorkspace);
372 reducedTriangle = qrDecomp.matrixQR().topRows(diagSize());
373 reducedTriangle.template triangularView<StrictlyLower>().setZero();
374 bid.compute(reducedTriangle);
376 bid.compute(copyWorkspace);
383 m_computed.topRows(diagSize()) = bid.bidiagonal().toDenseMatrix().transpose();
384 m_computed.template bottomRows<1>().setZero();
385 divide(0, diagSize() - 1, 0, 0, 0);
387 m_isInitialized =
true;
392 for (
int i = 0; i < diagSize(); i++) {
393 RealScalar a =
abs(m_computed.coeff(i, i));
394 m_singularValues.coeffRef(i) = a * scale;
395 if (a < considerZero) {
396 m_nonzeroSingularValues = i;
397 m_singularValues.tail(diagSize() - i - 1).setZero();
399 }
else if (i == diagSize() - 1) {
400 m_nonzeroSingularValues = i + 1;
407 copyUV(bid.householderV(), bid.householderU(), m_naiveV, m_naiveU);
409 copyUV(bid.householderU(), bid.householderV(), m_naiveU, m_naiveV);
412 if (m_isTranspose && computeV())
413 m_matrixV.applyOnTheLeft(qrDecomp.householderQ());
414 else if (!m_isTranspose && computeU())
415 m_matrixU.applyOnTheLeft(qrDecomp.householderQ());
418 m_isInitialized =
true;
422template <
typename MatrixType,
int Options>
423template <
typename HouseholderU,
typename HouseholderV,
typename NaiveU,
typename NaiveV>
424void BDCSVD<MatrixType, Options>::copyUV(
const HouseholderU& householderU,
const HouseholderV& householderV,
425 const NaiveU& naiveU,
const NaiveV& naiveV) {
428 Index Ucols = m_computeThinU ? diagSize() : rows();
429 m_matrixU = MatrixX::Identity(rows(), Ucols);
430 m_matrixU.topLeftCorner(diagSize(), diagSize()) =
431 naiveV.template cast<Scalar>().topLeftCorner(diagSize(), diagSize());
434 m_matrixU.topLeftCorner(householderU.cols(), diagSize()).applyOnTheLeft(householderU);
436 m_matrixU.applyOnTheLeft(householderU);
439 Index Vcols = m_computeThinV ? diagSize() : cols();
440 m_matrixV = MatrixX::Identity(cols(), Vcols);
441 m_matrixV.topLeftCorner(diagSize(), diagSize()) =
442 naiveU.template cast<Scalar>().topLeftCorner(diagSize(), diagSize());
445 m_matrixV.topLeftCorner(householderV.cols(), diagSize()).applyOnTheLeft(householderV);
447 m_matrixV.applyOnTheLeft(householderV);
459template <
typename MatrixType,
int Options>
470 Index k1 = 0, k2 = 0;
471 for (
Index j = 0; j < n; ++j) {
472 if ((A.col(j).head(n1).array() !=
Literal(0)).any()) {
473 A1.col(k1) = A.col(j).head(n1);
474 B1.row(k1) = B.row(j);
477 if ((A.col(j).tail(n2).array() !=
Literal(0)).any()) {
478 A2.col(k2) = A.col(j).tail(n2);
479 B2.row(k2) = B.row(j);
484 A.topRows(n1).noalias() = A1.leftCols(k1) * B1.topRows(k1);
485 A.bottomRows(n2).noalias() = A2.leftCols(k2) * B2.topRows(k2);
488 tmp.noalias() = A * B;
493template <
typename MatrixType,
int Options>
494template <
typename SVDType>
495void BDCSVD<MatrixType, Options>::computeBaseCase(SVDType& svd,
Index n,
Index firstCol,
Index firstRowW,
497 svd.compute(m_computed.block(firstCol, firstCol, n + 1, n));
501 m_naiveU.block(firstCol, firstCol, n + 1, n + 1).real() = svd.matrixU();
503 m_naiveU.row(0).segment(firstCol, n + 1).real() = svd.matrixU().row(0);
504 m_naiveU.row(1).segment(firstCol, n + 1).real() = svd.matrixU().row(n);
506 if (m_compV) m_naiveV.block(firstRowW, firstColW, n, n).real() = svd.matrixV();
507 m_computed.block(firstCol + shift, firstCol + shift, n + 1, n).setZero();
508 m_computed.diagonal().segment(firstCol + shift, n) = svd.singularValues().head(n);
524template <
typename MatrixType,
int Options>
530 const Index n = lastCol - firstCol + 1;
531 const Index k = n / 2;
532 const RealScalar considerZero = (std::numeric_limits<RealScalar>::min)();
536 RealScalar lambda, phi, c0, s0;
540 if (n < m_algoswap) {
544 computeBaseCase(baseSvd, n, firstCol, firstRowW, firstColW, shift);
547 computeBaseCase(baseSvd, n, firstCol, firstRowW, firstColW, shift);
552 alphaK = m_computed(firstCol + k, firstCol + k);
553 betaK = m_computed(firstCol + k + 1, firstCol + k);
557 divide(k + 1 + firstCol, lastCol, k + 1 + firstRowW, k + 1 + firstColW, shift);
559 divide(firstCol, k - 1 + firstCol, firstRowW, firstColW + 1, shift + 1);
563 lambda = m_naiveU(firstCol + k, firstCol + k);
564 phi = m_naiveU(firstCol + k + 1, lastCol + 1);
566 lambda = m_naiveU(1, firstCol + k);
567 phi = m_naiveU(0, lastCol + 1);
569 r0 =
sqrt((
abs(alphaK * lambda) *
abs(alphaK * lambda)) +
abs(betaK * phi) *
abs(betaK * phi));
571 l = m_naiveU.row(firstCol + k).segment(firstCol, k);
572 f = m_naiveU.row(firstCol + k + 1).segment(firstCol + k + 1, n - k - 1);
574 l = m_naiveU.row(1).segment(firstCol, k);
575 f = m_naiveU.row(0).segment(firstCol + k + 1, n - k - 1);
577 if (m_compV) m_naiveV(firstRowW + k, firstColW) =
Literal(1);
578 if (r0 < considerZero) {
582 c0 = alphaK * lambda / r0;
583 s0 = betaK * phi / r0;
586#ifdef EIGEN_BDCSVD_SANITY_CHECKS
587 eigen_internal_assert(m_naiveU.allFinite());
588 eigen_internal_assert(m_naiveV.allFinite());
589 eigen_internal_assert(m_computed.allFinite());
593 MatrixXr q1(m_naiveU.col(firstCol + k).segment(firstCol, k + 1));
595 for (
Index i = firstCol + k - 1; i >= firstCol; i--)
596 m_naiveU.col(i + 1).segment(firstCol, k + 1) = m_naiveU.col(i).segment(firstCol, k + 1);
598 m_naiveU.col(firstCol).segment(firstCol, k + 1) = (q1 * c0);
600 m_naiveU.col(lastCol + 1).segment(firstCol, k + 1) = (q1 * (-s0));
602 m_naiveU.col(firstCol).segment(firstCol + k + 1, n - k) =
603 m_naiveU.col(lastCol + 1).segment(firstCol + k + 1, n - k) * s0;
605 m_naiveU.col(lastCol + 1).segment(firstCol + k + 1, n - k) *= c0;
607 RealScalar q1 = m_naiveU(0, firstCol + k);
609 for (
Index i = firstCol + k - 1; i >= firstCol; i--) m_naiveU(0, i + 1) = m_naiveU(0, i);
611 m_naiveU(0, firstCol) = (q1 * c0);
613 m_naiveU(0, lastCol + 1) = (q1 * (-s0));
615 m_naiveU(1, firstCol) = m_naiveU(1, lastCol + 1) * s0;
617 m_naiveU(1, lastCol + 1) *= c0;
618 m_naiveU.row(1).segment(firstCol + 1, k).setZero();
619 m_naiveU.row(0).segment(firstCol + k + 1, n - k - 1).setZero();
622#ifdef EIGEN_BDCSVD_SANITY_CHECKS
623 eigen_internal_assert(m_naiveU.allFinite());
624 eigen_internal_assert(m_naiveV.allFinite());
625 eigen_internal_assert(m_computed.allFinite());
628 m_computed(firstCol + shift, firstCol + shift) = r0;
629 m_computed.col(firstCol + shift).segment(firstCol + shift + 1, k) = alphaK * l.transpose().real();
630 m_computed.col(firstCol + shift).segment(firstCol + shift + k + 1, n - k - 1) = betaK * f.transpose().real();
632#ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
633 ArrayXr tmp1 = (m_computed.block(firstCol + shift, firstCol + shift, n, n)).jacobiSvd().singularValues();
636 deflation(firstCol, lastCol, k, firstRowW, firstColW, shift);
637#ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
638 ArrayXr tmp2 = (m_computed.block(firstCol + shift, firstCol + shift, n, n)).jacobiSvd().singularValues();
639 std::cout <<
"\n\nj1 = " << tmp1.transpose().format(bdcsvdfmt) <<
"\n";
640 std::cout <<
"j2 = " << tmp2.transpose().format(bdcsvdfmt) <<
"\n\n";
641 std::cout <<
"err: " << ((tmp1 - tmp2).abs() > 1e-12 * tmp2.abs()).transpose() <<
"\n";
642 static int count = 0;
643 std::cout <<
"# " << ++count <<
"\n\n";
644 eigen_internal_assert((tmp1 - tmp2).matrix().norm() < 1e-14 * tmp2.matrix().norm());
650 MatrixXr UofSVD, VofSVD;
652 computeSVDofM(firstCol + shift, n, UofSVD, singVals, VofSVD);
654#ifdef EIGEN_BDCSVD_SANITY_CHECKS
655 eigen_internal_assert(UofSVD.allFinite());
656 eigen_internal_assert(VofSVD.allFinite());
660 structured_update(m_naiveU.block(firstCol, firstCol, n + 1, n + 1), UofSVD, (n + 2) / 2);
663 tmp.noalias() = m_naiveU.middleCols(firstCol, n + 1) * UofSVD;
664 m_naiveU.middleCols(firstCol, n + 1) = tmp;
667 if (m_compV) structured_update(m_naiveV.block(firstRowW, firstColW, n, n), VofSVD, (n + 1) / 2);
669#ifdef EIGEN_BDCSVD_SANITY_CHECKS
670 eigen_internal_assert(m_naiveU.allFinite());
671 eigen_internal_assert(m_naiveV.allFinite());
672 eigen_internal_assert(m_computed.allFinite());
675 m_computed.block(firstCol + shift, firstCol + shift, n, n).setZero();
676 m_computed.block(firstCol + shift, firstCol + shift, n, n).diagonal() = singVals;
688template <
typename MatrixType,
int Options>
689void BDCSVD<MatrixType, Options>::computeSVDofM(
Index firstCol,
Index n, MatrixXr& U, VectorType& singVals,
691 const RealScalar considerZero = (std::numeric_limits<RealScalar>::min)();
693 ArrayRef col0 = m_computed.col(firstCol).segment(firstCol, n);
694 m_workspace.head(n) = m_computed.block(firstCol, firstCol, n, n).diagonal();
695 ArrayRef diag = m_workspace.head(n);
700 U.resize(n + 1, n + 1);
701 if (m_compV) V.resize(n, n);
703#ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
704 if (col0.hasNaN() || diag.hasNaN()) std::cout <<
"\n\nHAS NAN\n\n";
711 while (actual_n > 1 && numext::is_exactly_zero(diag(actual_n - 1))) {
713 eigen_internal_assert(numext::is_exactly_zero(col0(actual_n)));
716 for (
Index k = 0; k < actual_n; ++k)
717 if (
abs(col0(k)) > considerZero) m_workspaceI(m++) = k;
724#ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
725 std::cout <<
"computeSVDofM using:\n";
726 std::cout <<
" z: " << col0.transpose() <<
"\n";
727 std::cout <<
" d: " << diag.transpose() <<
"\n";
731 computeSingVals(col0, diag, perm, singVals, shifts, mus);
733#ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
735 << (m_computed.block(firstCol, firstCol, n, n)).jacobiSvd().singularValues().transpose().reverse()
737 std::cout <<
" sing-val: " << singVals.transpose() <<
"\n";
738 std::cout <<
" mu: " << mus.transpose() <<
"\n";
739 std::cout <<
" shift: " << shifts.transpose() <<
"\n";
742 std::cout <<
"\n\n mus: " << mus.head(actual_n).transpose() <<
"\n\n";
743 std::cout <<
" check1 (expect0) : "
744 << ((singVals.array() - (shifts + mus)) / singVals.array()).head(actual_n).transpose() <<
"\n\n";
745 eigen_internal_assert((((singVals.array() - (shifts + mus)) / singVals.array()).head(actual_n) >= 0).all());
746 std::cout <<
" check2 (>0) : " << ((singVals.array() - diag) / singVals.array()).head(actual_n).transpose()
748 eigen_internal_assert((((singVals.array() - diag) / singVals.array()).head(actual_n) >= 0).all());
752#ifdef EIGEN_BDCSVD_SANITY_CHECKS
753 eigen_internal_assert(singVals.allFinite());
754 eigen_internal_assert(mus.allFinite());
755 eigen_internal_assert(shifts.allFinite());
759 perturbCol0(col0, diag, perm, singVals, shifts, mus, zhat);
760#ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
761 std::cout <<
" zhat: " << zhat.transpose() <<
"\n";
764#ifdef EIGEN_BDCSVD_SANITY_CHECKS
765 eigen_internal_assert(zhat.allFinite());
768 computeSingVecs(zhat, diag, perm, singVals, shifts, mus, U, V);
770#ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
771 std::cout <<
"U^T U: " << (U.transpose() * U - MatrixXr(MatrixXr::Identity(U.cols(), U.cols()))).norm() <<
"\n";
772 std::cout <<
"V^T V: " << (V.transpose() * V - MatrixXr(MatrixXr::Identity(V.cols(), V.cols()))).norm() <<
"\n";
775#ifdef EIGEN_BDCSVD_SANITY_CHECKS
776 eigen_internal_assert(m_naiveU.allFinite());
777 eigen_internal_assert(m_naiveV.allFinite());
778 eigen_internal_assert(m_computed.allFinite());
779 eigen_internal_assert(U.allFinite());
780 eigen_internal_assert(V.allFinite());
788 for (
Index i = 0; i < actual_n - 1; ++i) {
789 if (singVals(i) > singVals(i + 1)) {
791 swap(singVals(i), singVals(i + 1));
792 U.col(i).swap(U.col(i + 1));
793 if (m_compV) V.col(i).swap(V.col(i + 1));
797#ifdef EIGEN_BDCSVD_SANITY_CHECKS
799 bool singular_values_sorted =
800 (((singVals.segment(1, actual_n - 1) - singVals.head(actual_n - 1))).array() >= 0).all();
801 if (!singular_values_sorted)
802 std::cout <<
"Singular values are not sorted: " << singVals.segment(1, actual_n).transpose() <<
"\n";
803 eigen_internal_assert(singular_values_sorted);
809 singVals.head(actual_n).reverseInPlace();
810 U.leftCols(actual_n).rowwise().reverseInPlace();
811 if (m_compV) V.leftCols(actual_n).rowwise().reverseInPlace();
813#ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
815 std::cout <<
" * j: " << jsvd.singularValues().transpose() <<
"\n\n";
816 std::cout <<
" * sing-val: " << singVals.transpose() <<
"\n";
821template <
typename MatrixType,
int Options>
822typename BDCSVD<MatrixType, Options>::RealScalar BDCSVD<MatrixType, Options>::secularEq(
823 RealScalar mu,
const ArrayRef& col0,
const ArrayRef& diag,
const IndicesRef& perm,
const ArrayRef& diagShifted,
825 Index m = perm.size();
827 for (
Index i = 0; i < m; ++i) {
831 res += (col0(j) / (diagShifted(j) - mu)) * (col0(j) / (diag(j) + shift + mu));
836template <
typename MatrixType,
int Options>
837void BDCSVD<MatrixType, Options>::computeSingVals(
const ArrayRef& col0,
const ArrayRef& diag,
const IndicesRef& perm,
838 VectorType& singVals, ArrayRef shifts, ArrayRef mus) {
843 Index n = col0.size();
847 while (actual_n > 1 && numext::is_exactly_zero(col0(actual_n - 1))) --actual_n;
849 for (
Index k = 0; k < n; ++k) {
850 if (numext::is_exactly_zero(col0(k)) || actual_n == 1) {
853 singVals(k) = k == 0 ? col0(0) : diag(k);
855 shifts(k) = k == 0 ? col0(0) : diag(k);
860 RealScalar left = diag(k);
862 if (k == actual_n - 1)
863 right = (diag(actual_n - 1) + col0.matrix().norm());
869 while (numext::is_exactly_zero(col0(l))) {
871 eigen_internal_assert(l < actual_n);
877 RealScalar mid = left + (right - left) /
Literal(2);
878 RealScalar fMid = secularEq(mid, col0, diag, perm, diag,
Literal(0));
879#ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
880 std::cout <<
"right-left = " << right - left <<
"\n";
884 std::cout <<
" = " << secularEq(left + RealScalar(0.000001) * (right - left), col0, diag, perm, diag, 0) <<
" "
885 << secularEq(left + RealScalar(0.1) * (right - left), col0, diag, perm, diag, 0) <<
" "
886 << secularEq(left + RealScalar(0.2) * (right - left), col0, diag, perm, diag, 0) <<
" "
887 << secularEq(left + RealScalar(0.3) * (right - left), col0, diag, perm, diag, 0) <<
" "
888 << secularEq(left + RealScalar(0.4) * (right - left), col0, diag, perm, diag, 0) <<
" "
889 << secularEq(left + RealScalar(0.49) * (right - left), col0, diag, perm, diag, 0) <<
" "
890 << secularEq(left + RealScalar(0.5) * (right - left), col0, diag, perm, diag, 0) <<
" "
891 << secularEq(left + RealScalar(0.51) * (right - left), col0, diag, perm, diag, 0) <<
" "
892 << secularEq(left + RealScalar(0.6) * (right - left), col0, diag, perm, diag, 0) <<
" "
893 << secularEq(left + RealScalar(0.7) * (right - left), col0, diag, perm, diag, 0) <<
" "
894 << secularEq(left + RealScalar(0.8) * (right - left), col0, diag, perm, diag, 0) <<
" "
895 << secularEq(left + RealScalar(0.9) * (right - left), col0, diag, perm, diag, 0) <<
" "
896 << secularEq(left + RealScalar(0.999999) * (right - left), col0, diag, perm, diag, 0) <<
"\n";
898 RealScalar shift = (k == actual_n - 1 || fMid >
Literal(0)) ? left : right;
901 Map<ArrayXr> diagShifted(m_workspace.data() + 4 * n, n);
902 diagShifted = diag - shift;
904 if (k != actual_n - 1) {
906 RealScalar midShifted = (right - left) / RealScalar(2);
908 if (numext::equal_strict(shift, right)) midShifted = -midShifted;
909 RealScalar fMidShifted = secularEq(midShifted, col0, diag, perm, diagShifted, shift);
910 if (fMidShifted > 0) {
912 shift = fMidShifted >
Literal(0) ? left : right;
913 diagShifted = diag - shift;
918 RealScalar muPrev, muCur;
920 if (numext::equal_strict(shift, left)) {
921 muPrev = (right - left) * RealScalar(0.1);
922 if (k == actual_n - 1)
923 muCur = right - left;
925 muCur = (right - left) * RealScalar(0.5);
927 muPrev = -(right - left) * RealScalar(0.1);
928 muCur = -(right - left) * RealScalar(0.5);
931 RealScalar fPrev = secularEq(muPrev, col0, diag, perm, diagShifted, shift);
932 RealScalar fCur = secularEq(muCur, col0, diag, perm, diagShifted, shift);
933 if (
abs(fPrev) <
abs(fCur)) {
940 bool useBisection = fPrev * fCur >
Literal(0);
941 while (!numext::is_exactly_zero(fCur) &&
942 abs(muCur - muPrev) >
943 Literal(8) * NumTraits<RealScalar>::epsilon() * numext::maxi<RealScalar>(
abs(muCur),
abs(muPrev)) &&
944 abs(fCur - fPrev) > NumTraits<RealScalar>::epsilon() && !useBisection) {
948 RealScalar a = (fCur - fPrev) / (
Literal(1) / muCur -
Literal(1) / muPrev);
949 RealScalar b = fCur - a / muCur;
951 RealScalar muZero = -a / b;
952 RealScalar fZero = secularEq(muZero, col0, diag, perm, diagShifted, shift);
954#ifdef EIGEN_BDCSVD_SANITY_CHECKS
955 eigen_internal_assert((numext::isfinite)(fZero));
964 if (numext::equal_strict(shift, left) && (muCur <
Literal(0) || muCur > right - left)) useBisection =
true;
965 if (numext::equal_strict(shift, right) && (muCur < -(right - left) || muCur >
Literal(0))) useBisection =
true;
966 if (
abs(fCur) >
abs(fPrev)) useBisection =
true;
971#ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
972 std::cout <<
"useBisection for k = " << k <<
", actual_n = " << actual_n <<
"\n";
974 RealScalar leftShifted, rightShifted;
976 if (numext::equal_strict(shift, left)) {
980 numext::maxi<RealScalar>((std::numeric_limits<RealScalar>::min)(),
981 Literal(2) *
abs(col0(k)) /
sqrt((std::numeric_limits<RealScalar>::max)()));
984 eigen_internal_assert(
985 (numext::isfinite)((col0(k) / leftShifted) * (col0(k) / (diag(k) + shift + leftShifted))));
988 rightShifted = (k == actual_n - 1)
990 : ((right - left) * RealScalar(0.51));
992 leftShifted = -(right - left) * RealScalar(0.51);
994 rightShifted = -numext::maxi<RealScalar>((std::numeric_limits<RealScalar>::min)(),
995 abs(col0(k + 1)) /
sqrt((std::numeric_limits<RealScalar>::max)()));
997 rightShifted = -(std::numeric_limits<RealScalar>::min)();
1000 RealScalar fLeft = secularEq(leftShifted, col0, diag, perm, diagShifted, shift);
1001 eigen_internal_assert(fLeft <
Literal(0));
1003#if defined EIGEN_BDCSVD_DEBUG_VERBOSE || defined EIGEN_BDCSVD_SANITY_CHECKS || defined EIGEN_INTERNAL_DEBUGGING
1004 RealScalar fRight = secularEq(rightShifted, col0, diag, perm, diagShifted, shift);
1007#ifdef EIGEN_BDCSVD_SANITY_CHECKS
1008 if (!(numext::isfinite)(fLeft))
1009 std::cout <<
"f(" << leftShifted <<
") =" << fLeft <<
" ; " << left <<
" " << shift <<
" " << right <<
"\n";
1010 eigen_internal_assert((numext::isfinite)(fLeft));
1012 if (!(numext::isfinite)(fRight))
1013 std::cout <<
"f(" << rightShifted <<
") =" << fRight <<
" ; " << left <<
" " << shift <<
" " << right <<
"\n";
1017#ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
1018 if (!(fLeft * fRight < 0)) {
1019 std::cout <<
"f(leftShifted) using leftShifted=" << leftShifted
1020 <<
" ; diagShifted(1:10):" << diagShifted.head(10).transpose() <<
"\n ; "
1021 <<
"left==shift=" << bool(left == shift) <<
" ; left-shift = " << (left - shift) <<
"\n";
1022 std::cout <<
"k=" << k <<
", " << fLeft <<
" * " << fRight <<
" == " << fLeft * fRight <<
" ; "
1023 <<
"[" << left <<
" .. " << right <<
"] -> [" << leftShifted <<
" " << rightShifted
1024 <<
"], shift=" << shift <<
" , f(right)=" << secularEq(0, col0, diag, perm, diagShifted, shift)
1025 <<
" == " << secularEq(right, col0, diag, perm, diag, 0) <<
" == " << fRight <<
"\n";
1028 eigen_internal_assert(fLeft * fRight <
Literal(0));
1031 while (rightShifted - leftShifted >
Literal(2) * NumTraits<RealScalar>::epsilon() *
1032 numext::maxi<RealScalar>(
abs(leftShifted),
abs(rightShifted))) {
1033 RealScalar midShifted = (leftShifted + rightShifted) /
Literal(2);
1034 fMid = secularEq(midShifted, col0, diag, perm, diagShifted, shift);
1035 eigen_internal_assert((numext::isfinite)(fMid));
1037 if (fLeft * fMid <
Literal(0)) {
1038 rightShifted = midShifted;
1040 leftShifted = midShifted;
1044 muCur = (leftShifted + rightShifted) /
Literal(2);
1050 muCur = (right - left) * RealScalar(0.5);
1052 if (numext::equal_strict(shift, right)) muCur = -muCur;
1056 singVals[k] = shift + muCur;
1060#ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
1062 std::cout <<
"found " << singVals[k] <<
" == " << shift <<
" + " << muCur <<
" from " << diag(k) <<
" .. "
1063 << diag(k + 1) <<
"\n";
1065#ifdef EIGEN_BDCSVD_SANITY_CHECKS
1066 eigen_internal_assert(k == 0 || singVals[k] >= singVals[k - 1]);
1067 eigen_internal_assert(singVals[k] >= diag(k));
1079template <
typename MatrixType,
int Options>
1080void BDCSVD<MatrixType, Options>::perturbCol0(
const ArrayRef& col0,
const ArrayRef& diag,
const IndicesRef& perm,
1081 const VectorType& singVals,
const ArrayRef& shifts,
const ArrayRef& mus,
1084 Index n = col0.size();
1085 Index m = perm.size();
1090 Index lastIdx = perm(m - 1);
1092 for (
Index k = 0; k < n; ++k) {
1093 if (numext::is_exactly_zero(col0(k)))
1097 RealScalar dk = diag(k);
1098 RealScalar prod = (singVals(lastIdx) + dk) * (mus(lastIdx) + (shifts(lastIdx) - dk));
1099#ifdef EIGEN_BDCSVD_SANITY_CHECKS
1101 std::cout <<
"k = " << k <<
" ; z(k)=" << col0(k) <<
", diag(k)=" << dk <<
"\n";
1102 std::cout <<
"prod = "
1103 <<
"(" << singVals(lastIdx) <<
" + " << dk <<
") * (" << mus(lastIdx) <<
" + (" << shifts(lastIdx)
1104 <<
" - " << dk <<
"))"
1106 std::cout <<
" = " << singVals(lastIdx) + dk <<
" * " << mus(lastIdx) + (shifts(lastIdx) - dk) <<
"\n";
1108 eigen_internal_assert(prod >= 0);
1111 for (
Index l = 0; l < m; ++l) {
1114#ifdef EIGEN_BDCSVD_SANITY_CHECKS
1115 if (i >= k && (l == 0 || l - 1 >= m)) {
1116 std::cout <<
"Error in perturbCol0\n";
1117 std::cout <<
" " << k <<
"/" << n <<
" " << l <<
"/" << m <<
" " << i <<
"/" << n <<
" ; " << col0(k)
1118 <<
" " << diag(k) <<
" "
1120 std::cout <<
" " << diag(i) <<
"\n";
1121 Index j = (i < k ) ? i : perm(l - 1);
1123 <<
"j=" << j <<
"\n";
1126 Index j = i < k ? i : l > 0 ? perm(l - 1) : i;
1127#ifdef EIGEN_BDCSVD_SANITY_CHECKS
1129 std::cout <<
"k=" << k <<
", i=" << i <<
", l=" << l <<
", perm.size()=" << perm.size() <<
"\n";
1131 eigen_internal_assert(dk !=
Literal(0) || diag(i) !=
Literal(0));
1133 prod *= ((singVals(j) + dk) / ((diag(i) + dk))) * ((mus(j) + (shifts(j) - dk)) / ((diag(i) - dk)));
1134#ifdef EIGEN_BDCSVD_SANITY_CHECKS
1135 eigen_internal_assert(prod >= 0);
1137#ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
1139 numext::abs(((singVals(j) + dk) * (mus(j) + (shifts(j) - dk))) / ((diag(i) + dk) * (diag(i) - dk)) - 1) >
1142 << ((singVals(j) + dk) * (mus(j) + (shifts(j) - dk))) / ((diag(i) + dk) * (diag(i) - dk))
1143 <<
" == (" << (singVals(j) + dk) <<
" * " << (mus(j) + (shifts(j) - dk)) <<
") / ("
1144 << (diag(i) + dk) <<
" * " << (diag(i) - dk) <<
")\n";
1148#ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
1149 std::cout <<
"zhat(" << k <<
") = sqrt( " << prod <<
") ; " << (singVals(lastIdx) + dk) <<
" * "
1150 << mus(lastIdx) + shifts(lastIdx) <<
" - " << dk <<
"\n";
1152 RealScalar tmp =
sqrt(prod);
1153#ifdef EIGEN_BDCSVD_SANITY_CHECKS
1154 eigen_internal_assert((numext::isfinite)(tmp));
1156 zhat(k) = col0(k) >
Literal(0) ? RealScalar(tmp) : RealScalar(-tmp);
1162template <
typename MatrixType,
int Options>
1163void BDCSVD<MatrixType, Options>::computeSingVecs(
const ArrayRef& zhat,
const ArrayRef& diag,
const IndicesRef& perm,
1164 const VectorType& singVals,
const ArrayRef& shifts,
1165 const ArrayRef& mus, MatrixXr& U, MatrixXr& V) {
1166 Index n = zhat.size();
1167 Index m = perm.size();
1169 for (
Index k = 0; k < n; ++k) {
1170 if (numext::is_exactly_zero(zhat(k))) {
1171 U.col(k) = VectorType::Unit(n + 1, k);
1172 if (m_compV) V.col(k) = VectorType::Unit(n, k);
1175 for (
Index l = 0; l < m; ++l) {
1177 U(i, k) = zhat(i) / (((diag(i) - shifts(k)) - mus(k))) / ((diag(i) + singVals[k]));
1180 U.col(k).normalize();
1184 for (
Index l = 1; l < m; ++l) {
1186 V(i, k) = diag(i) * zhat(i) / (((diag(i) - shifts(k)) - mus(k))) / ((diag(i) + singVals[k]));
1189 V.col(k).normalize();
1193 U.col(n) = VectorType::Unit(n + 1, n);
1199template <
typename MatrixType,
int Options>
1200void BDCSVD<MatrixType, Options>::deflation43(
Index firstCol,
Index shift,
Index i,
Index size) {
1204 Index start = firstCol + shift;
1205 RealScalar c = m_computed(start, start);
1206 RealScalar s = m_computed(start + i, start);
1207 RealScalar r = numext::hypot(c, s);
1208 if (numext::is_exactly_zero(r)) {
1209 m_computed(start + i, start + i) =
Literal(0);
1212 m_computed(start, start) = r;
1213 m_computed(start + i, start) =
Literal(0);
1214 m_computed(start + i, start + i) =
Literal(0);
1218 m_naiveU.middleRows(firstCol, size + 1).applyOnTheRight(firstCol, firstCol + i, J);
1220 m_naiveU.applyOnTheRight(firstCol, firstCol + i, J);
1226template <
typename MatrixType,
int Options>
1227void BDCSVD<MatrixType, Options>::deflation44(
Index firstColu,
Index firstColm,
Index firstRowW,
Index firstColW,
1234 RealScalar s = m_computed(firstColm + i, firstColm);
1235 RealScalar c = m_computed(firstColm + j, firstColm);
1236 RealScalar r = numext::hypot(c, s);
1237#ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
1238 std::cout <<
"deflation 4.4: " << i <<
"," << j <<
" -> " << c <<
" " << s <<
" " << r <<
" ; "
1239 << m_computed(firstColm + i - 1, firstColm) <<
" " << m_computed(firstColm + i, firstColm) <<
" "
1240 << m_computed(firstColm + i + 1, firstColm) <<
" " << m_computed(firstColm + i + 2, firstColm) <<
"\n";
1241 std::cout << m_computed(firstColm + i - 1, firstColm + i - 1) <<
" " << m_computed(firstColm + i, firstColm + i)
1242 <<
" " << m_computed(firstColm + i + 1, firstColm + i + 1) <<
" "
1243 << m_computed(firstColm + i + 2, firstColm + i + 2) <<
"\n";
1245 if (numext::is_exactly_zero(r)) {
1246 m_computed(firstColm + j, firstColm + j) = m_computed(firstColm + i, firstColm + i);
1251 m_computed(firstColm + j, firstColm) = r;
1252 m_computed(firstColm + j, firstColm + j) = m_computed(firstColm + i, firstColm + i);
1253 m_computed(firstColm + i, firstColm) =
Literal(0);
1257 m_naiveU.middleRows(firstColu, size + 1).applyOnTheRight(firstColu + j, firstColu + i, J);
1259 m_naiveU.applyOnTheRight(firstColu + j, firstColu + i, J);
1260 if (m_compV) m_naiveV.middleRows(firstRowW, size).applyOnTheRight(firstColW + j, firstColW + i, J);
1264template <
typename MatrixType,
int Options>
1269 const Index length = lastCol + 1 - firstCol;
1275 const RealScalar considerZero = (std::numeric_limits<RealScalar>::min)();
1276 RealScalar maxDiag = diag.tail((std::max)(
Index(1), length - 1)).cwiseAbs().maxCoeff();
1277 RealScalar epsilon_strict = numext::maxi<RealScalar>(considerZero, NumTraits<RealScalar>::epsilon() * maxDiag);
1278 RealScalar epsilon_coarse =
1279 Literal(8) * NumTraits<RealScalar>::epsilon() * numext::maxi<RealScalar>(col0.cwiseAbs().maxCoeff(), maxDiag);
1281#ifdef EIGEN_BDCSVD_SANITY_CHECKS
1282 eigen_internal_assert(m_naiveU.allFinite());
1283 eigen_internal_assert(m_naiveV.allFinite());
1284 eigen_internal_assert(m_computed.allFinite());
1287#ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
1288 std::cout <<
"\ndeflate:" << diag.head(k + 1).transpose() <<
" | "
1289 << diag.segment(k + 1, length - k - 1).transpose() <<
"\n";
1293 if (diag(0) < epsilon_coarse) {
1294#ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
1295 std::cout <<
"deflation 4.1, because " << diag(0) <<
" < " << epsilon_coarse <<
"\n";
1297 diag(0) = epsilon_coarse;
1301 for (
Index i = 1; i < length; ++i)
1302 if (
abs(col0(i)) < epsilon_strict) {
1303#ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
1304 std::cout <<
"deflation 4.2, set z(" << i <<
") to zero because " <<
abs(col0(i)) <<
" < " << epsilon_strict
1305 <<
" (diag(" << i <<
")=" << diag(i) <<
")\n";
1311 for (
Index i = 1; i < length; i++)
1312 if (diag(i) < epsilon_coarse) {
1313#ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
1314 std::cout <<
"deflation 4.3, cancel z(" << i <<
")=" << col0(i) <<
" because diag(" << i <<
")=" << diag(i)
1315 <<
" < " << epsilon_coarse <<
"\n";
1317 deflation43(firstCol, shift, i, length);
1320#ifdef EIGEN_BDCSVD_SANITY_CHECKS
1321 eigen_internal_assert(m_naiveU.allFinite());
1322 eigen_internal_assert(m_naiveV.allFinite());
1323 eigen_internal_assert(m_computed.allFinite());
1325#ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
1326 std::cout <<
"to be sorted: " << diag.transpose() <<
"\n\n";
1327 std::cout <<
" : " << col0.transpose() <<
"\n\n";
1332 const bool total_deflation = (col0.tail(length - 1).array().
abs() < considerZero).all();
1336 Index* permutation = m_workspaceI.data();
1342 for (
Index i = 1; i < length; ++i)
1343 if (diag(i) < considerZero) permutation[p++] = i;
1345 Index i = 1, j = k + 1;
1346 for (; p < length; ++p) {
1348 permutation[p] = j++;
1349 else if (j >= length)
1350 permutation[p] = i++;
1351 else if (diag(i) < diag(j))
1352 permutation[p] = j++;
1354 permutation[p] = i++;
1359 if (total_deflation) {
1360 for (
Index i = 1; i < length; ++i) {
1361 Index pi = permutation[i];
1362 if (diag(pi) < considerZero || diag(0) < diag(pi))
1363 permutation[i - 1] = permutation[i];
1365 permutation[i - 1] = 0;
1372 Index* realInd = m_workspaceI.data() + length;
1373 Index* realCol = m_workspaceI.data() + 2 * length;
1375 for (
int pos = 0; pos < length; pos++) {
1380 for (
Index i = total_deflation ? 0 : 1; i < length; i++) {
1381 const Index pi = permutation[length - (total_deflation ? i + 1 : i)];
1382 const Index J = realCol[pi];
1386 swap(diag(i), diag(J));
1387 if (i != 0 && J != 0)
swap(col0(i), col0(J));
1391 m_naiveU.col(firstCol + i)
1392 .segment(firstCol, length + 1)
1393 .swap(m_naiveU.col(firstCol + J).segment(firstCol, length + 1));
1395 m_naiveU.col(firstCol + i).segment(0, 2).swap(m_naiveU.col(firstCol + J).segment(0, 2));
1397 m_naiveV.col(firstColW + i)
1398 .segment(firstRowW, length)
1399 .swap(m_naiveV.col(firstColW + J).segment(firstRowW, length));
1402 const Index realI = realInd[i];
1409#ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
1410 std::cout <<
"sorted: " << diag.transpose().format(bdcsvdfmt) <<
"\n";
1411 std::cout <<
" : " << col0.transpose() <<
"\n\n";
1416 Index i = length - 1;
1418 while (i > 0 && (diag(i) < considerZero ||
abs(col0(i)) < considerZero)) --i;
1421 if ((diag(i) - diag(i - 1)) < epsilon_strict) {
1422#ifdef EIGEN_BDCSVD_DEBUG_VERBOSE
1423 std::cout <<
"deflation 4.4 with i = " << i <<
" because " << diag(i) <<
" - " << diag(i - 1)
1424 <<
" == " << (diag(i) - diag(i - 1)) <<
" < " << epsilon_strict <<
"\n";
1426 eigen_internal_assert(
abs(diag(i) - diag(i - 1)) < epsilon_coarse &&
1427 " diagonal entries are not properly sorted");
1428 deflation44(firstCol, firstCol + shift, firstRowW, firstColW, i, i - 1, length);
1432#ifdef EIGEN_BDCSVD_SANITY_CHECKS
1433 for (
Index j = 2; j < length; ++j) eigen_internal_assert(diag(j - 1) <= diag(j) ||
abs(diag(j)) < considerZero);
1436#ifdef EIGEN_BDCSVD_SANITY_CHECKS
1437 eigen_internal_assert(m_naiveU.allFinite());
1438 eigen_internal_assert(m_naiveV.allFinite());
1439 eigen_internal_assert(m_computed.allFinite());
1449template <
typename Derived>
1450template <
int Options>
1461template <
typename Derived>
1462template <
int Options>
1464 unsigned int computationOptions)
const {
General-purpose arrays with easy API for coefficient-wise operations.
Definition Array.h:48
class Bidiagonal Divide and Conquer SVD
Definition BDCSVD.h:85
BDCSVD()
Default Constructor.
Definition BDCSVD.h:130
BDCSVD(const MatrixBase< Derived > &matrix)
Constructor performing the decomposition of given matrix, using the custom options specified with the...
Definition BDCSVD.h:170
BDCSVD & compute(const MatrixBase< Derived > &matrix)
Method performing the decomposition of given matrix. Computes Thin/Full unitaries U/V if specified us...
Definition BDCSVD.h:201
BDCSVD(Index rows, Index cols)
Default Constructor with memory preallocation.
Definition BDCSVD.h:138
Expression of a fixed-size or dynamic-size block.
Definition Block.h:110
Expression of a diagonal/subdiagonal/superdiagonal in a matrix.
Definition Diagonal.h:68
Householder QR decomposition of a matrix.
Definition HouseholderQR.h:59
Rotation given by a cosine-sine pair.
Definition Jacobi.h:38
Two-sided Jacobi SVD decomposition of a rectangular matrix.
Definition JacobiSVD.h:500
A matrix or vector expression mapping an existing array of data.
Definition Map.h:96
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
bool computeV() const
Definition SVDBase.h:275
bool computeU() const
Definition SVDBase.h:273
SVDBase()
Definition SVDBase.h:349
Expression of a fixed-size or dynamic-size sub-vector.
Definition VectorBlock.h:58
@ Aligned
Definition Constants.h:242
@ DisableQRDecomposition
Definition Constants.h:429
@ InvalidInput
Definition Constants.h:447
@ Success
Definition Constants.h:440
@ NoConvergence
Definition Constants.h:444
Matrix< Type, Dynamic, Dynamic > MatrixX
[c++11] Dynamic×Dynamic matrix of type Type.
Definition Matrix.h:514
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)
EIGEN_DEPRECATED_WITH_REASON("Initialization is no longer needed.") inline void initParallel()
Definition Parallelizer.h:50
std::enable_if_t< std::is_base_of< DenseBase< std::decay_t< DerivedA > >, std::decay_t< DerivedA > >::value &&std::is_base_of< DenseBase< std::decay_t< DerivedB > >, std::decay_t< DerivedB > >::value, void > swap(DerivedA &&a, DerivedB &&b)
Definition DenseBase.h:667
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:82
const int Dynamic
Definition Constants.h:25
Eigen::Index Index
The interface type of indices.
Definition EigenBase.h:43
constexpr Index size() const noexcept
Definition EigenBase.h:64