11#ifndef EIGEN_JACOBISVD_H
12#define EIGEN_JACOBISVD_H
15#include "./InternalHeaderCheck.h"
23template <typename MatrixType, int Options, bool IsComplex = NumTraits<typename MatrixType::Scalar>::IsComplex>
24struct svd_precondition_2x2_block_to_be_real {};
33enum { PreconditionIfMoreColsThanRows, PreconditionIfMoreRowsThanCols };
35template <
typename MatrixType,
int QRPreconditioner,
int Case>
36struct qr_preconditioner_should_do_anything {
38 a = MatrixType::RowsAtCompileTime !=
Dynamic && MatrixType::ColsAtCompileTime !=
Dynamic &&
39 MatrixType::ColsAtCompileTime <= MatrixType::RowsAtCompileTime,
40 b = MatrixType::RowsAtCompileTime !=
Dynamic && MatrixType::ColsAtCompileTime !=
Dynamic &&
41 MatrixType::RowsAtCompileTime <= MatrixType::ColsAtCompileTime,
42 ret = !((QRPreconditioner ==
NoQRPreconditioner) || (Case == PreconditionIfMoreColsThanRows && bool(a)) ||
43 (Case == PreconditionIfMoreRowsThanCols && bool(b)))
47template <
typename MatrixType,
int Options,
int QRPreconditioner,
int Case,
48 bool DoAnything = qr_preconditioner_should_do_anything<MatrixType, QRPreconditioner, Case>::ret>
49struct qr_preconditioner_impl {};
51template <
typename MatrixType,
int Options,
int QRPreconditioner,
int Case>
52class qr_preconditioner_impl<MatrixType, Options, QRPreconditioner, Case, false> {
54 void allocate(
const JacobiSVD<MatrixType, Options>&) {}
55 template <
typename Xpr>
56 bool run(JacobiSVD<MatrixType, Options>&,
const Xpr&) {
63template <
typename MatrixType,
int Options>
67 typedef typename MatrixType::Scalar Scalar;
68 typedef JacobiSVD<MatrixType, Options> SVDType;
70 enum { WorkspaceSize = MatrixType::RowsAtCompileTime, MaxWorkspaceSize = MatrixType::MaxRowsAtCompileTime };
72 typedef Matrix<Scalar, 1, WorkspaceSize, RowMajor, 1, MaxWorkspaceSize> WorkspaceType;
74 void allocate(
const SVDType& svd) {
75 if (svd.rows() != m_qr.rows() || svd.cols() != m_qr.cols()) {
76 internal::destroy_at(&m_qr);
77 internal::construct_at(&m_qr, svd.rows(), svd.cols());
79 if (svd.m_computeFullU) m_workspace.resize(svd.rows());
81 template <
typename Xpr>
82 bool run(SVDType& svd,
const Xpr& matrix) {
83 if (matrix.rows() > matrix.cols()) {
85 svd.m_workMatrix = m_qr.matrixQR().block(0, 0, matrix.cols(), matrix.cols()).template triangularView<Upper>();
86 if (svd.m_computeFullU) m_qr.matrixQ().evalTo(svd.m_matrixU, m_workspace);
87 if (svd.computeV()) svd.m_matrixV = m_qr.colsPermutation();
94 typedef FullPivHouseholderQR<MatrixType> QRType;
96 WorkspaceType m_workspace;
99template <
typename MatrixType,
int Options>
103 typedef typename MatrixType::Scalar Scalar;
104 typedef JacobiSVD<MatrixType, Options> SVDType;
107 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
108 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
109 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
110 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
111 MatrixOptions = traits<MatrixType>::Options
114 typedef typename internal::make_proper_matrix_type<Scalar, ColsAtCompileTime, RowsAtCompileTime, MatrixOptions,
115 MaxColsAtCompileTime, MaxRowsAtCompileTime>::type
116 TransposeTypeWithSameStorageOrder;
118 void allocate(
const SVDType& svd) {
119 if (svd.cols() != m_qr.rows() || svd.rows() != m_qr.cols()) {
120 internal::destroy_at(&m_qr);
121 internal::construct_at(&m_qr, svd.cols(), svd.rows());
123 if (svd.m_computeFullV) m_workspace.resize(svd.cols());
125 template <
typename Xpr>
126 bool run(SVDType& svd,
const Xpr& matrix) {
127 if (matrix.cols() > matrix.rows()) {
128 m_qr.compute(matrix.adjoint());
130 m_qr.matrixQR().block(0, 0, matrix.rows(), matrix.rows()).template triangularView<Upper>().adjoint();
131 if (svd.m_computeFullV) m_qr.matrixQ().evalTo(svd.m_matrixV, m_workspace);
132 if (svd.computeU()) svd.m_matrixU = m_qr.colsPermutation();
139 typedef FullPivHouseholderQR<TransposeTypeWithSameStorageOrder> QRType;
141 typename plain_row_type<MatrixType>::type m_workspace;
146template <
typename MatrixType,
int Options>
150 typedef typename MatrixType::Scalar Scalar;
151 typedef JacobiSVD<MatrixType, Options> SVDType;
154 WorkspaceSize = internal::traits<SVDType>::MatrixUColsAtCompileTime,
155 MaxWorkspaceSize = internal::traits<SVDType>::MatrixUMaxColsAtCompileTime
158 typedef Matrix<Scalar, 1, WorkspaceSize, RowMajor, 1, MaxWorkspaceSize> WorkspaceType;
160 void allocate(
const SVDType& svd) {
161 if (svd.rows() != m_qr.rows() || svd.cols() != m_qr.cols()) {
162 internal::destroy_at(&m_qr);
163 internal::construct_at(&m_qr, svd.rows(), svd.cols());
165 if (svd.m_computeFullU)
166 m_workspace.resize(svd.rows());
167 else if (svd.m_computeThinU)
168 m_workspace.resize(svd.cols());
170 template <
typename Xpr>
171 bool run(SVDType& svd,
const Xpr& matrix) {
172 if (matrix.rows() > matrix.cols()) {
173 m_qr.compute(matrix);
174 svd.m_workMatrix = m_qr.matrixQR().block(0, 0, matrix.cols(), matrix.cols()).template triangularView<Upper>();
175 if (svd.m_computeFullU)
176 m_qr.householderQ().evalTo(svd.m_matrixU, m_workspace);
177 else if (svd.m_computeThinU) {
178 svd.m_matrixU.setIdentity(matrix.rows(), matrix.cols());
179 m_qr.householderQ().applyThisOnTheLeft(svd.m_matrixU, m_workspace);
181 if (svd.computeV()) svd.m_matrixV = m_qr.colsPermutation();
188 typedef ColPivHouseholderQR<MatrixType> QRType;
190 WorkspaceType m_workspace;
193template <
typename MatrixType,
int Options>
197 typedef typename MatrixType::Scalar Scalar;
198 typedef JacobiSVD<MatrixType, Options> SVDType;
201 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
202 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
203 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
204 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
205 MatrixOptions = internal::traits<MatrixType>::Options,
206 WorkspaceSize = internal::traits<SVDType>::MatrixVColsAtCompileTime,
207 MaxWorkspaceSize = internal::traits<SVDType>::MatrixVMaxColsAtCompileTime
210 typedef Matrix<Scalar, WorkspaceSize, 1, ColMajor, MaxWorkspaceSize, 1> WorkspaceType;
212 typedef typename internal::make_proper_matrix_type<Scalar, ColsAtCompileTime, RowsAtCompileTime, MatrixOptions,
213 MaxColsAtCompileTime, MaxRowsAtCompileTime>::type
214 TransposeTypeWithSameStorageOrder;
216 void allocate(
const SVDType& svd) {
217 if (svd.cols() != m_qr.rows() || svd.rows() != m_qr.cols()) {
218 internal::destroy_at(&m_qr);
219 internal::construct_at(&m_qr, svd.cols(), svd.rows());
221 if (svd.m_computeFullV)
222 m_workspace.resize(svd.cols());
223 else if (svd.m_computeThinV)
224 m_workspace.resize(svd.rows());
226 template <
typename Xpr>
227 bool run(SVDType& svd,
const Xpr& matrix) {
228 if (matrix.cols() > matrix.rows()) {
229 m_qr.compute(matrix.adjoint());
232 m_qr.matrixQR().block(0, 0, matrix.rows(), matrix.rows()).template triangularView<Upper>().adjoint();
233 if (svd.m_computeFullV)
234 m_qr.householderQ().evalTo(svd.m_matrixV, m_workspace);
235 else if (svd.m_computeThinV) {
236 svd.m_matrixV.setIdentity(matrix.cols(), matrix.rows());
237 m_qr.householderQ().applyThisOnTheLeft(svd.m_matrixV, m_workspace);
239 if (svd.computeU()) svd.m_matrixU = m_qr.colsPermutation();
246 typedef ColPivHouseholderQR<TransposeTypeWithSameStorageOrder> QRType;
248 WorkspaceType m_workspace;
253template <
typename MatrixType,
int Options>
256 typedef typename MatrixType::Scalar Scalar;
257 typedef JacobiSVD<MatrixType, Options> SVDType;
260 WorkspaceSize = internal::traits<SVDType>::MatrixUColsAtCompileTime,
261 MaxWorkspaceSize = internal::traits<SVDType>::MatrixUMaxColsAtCompileTime
264 typedef Matrix<Scalar, 1, WorkspaceSize, RowMajor, 1, MaxWorkspaceSize> WorkspaceType;
266 void allocate(
const SVDType& svd) {
267 if (svd.rows() != m_qr.rows() || svd.cols() != m_qr.cols()) {
268 internal::destroy_at(&m_qr);
269 internal::construct_at(&m_qr, svd.rows(), svd.cols());
271 if (svd.m_computeFullU)
272 m_workspace.resize(svd.rows());
273 else if (svd.m_computeThinU)
274 m_workspace.resize(svd.cols());
276 template <
typename Xpr>
277 bool run(SVDType& svd,
const Xpr& matrix) {
278 if (matrix.rows() > matrix.cols()) {
279 m_qr.compute(matrix);
280 svd.m_workMatrix = m_qr.matrixQR().block(0, 0, matrix.cols(), matrix.cols()).template triangularView<Upper>();
281 if (svd.m_computeFullU)
282 m_qr.householderQ().evalTo(svd.m_matrixU, m_workspace);
283 else if (svd.m_computeThinU) {
284 svd.m_matrixU.setIdentity(matrix.rows(), matrix.cols());
285 m_qr.householderQ().applyThisOnTheLeft(svd.m_matrixU, m_workspace);
287 if (svd.computeV()) svd.m_matrixV.setIdentity(matrix.cols(), matrix.cols());
294 typedef HouseholderQR<MatrixType> QRType;
296 WorkspaceType m_workspace;
299template <
typename MatrixType,
int Options>
302 typedef typename MatrixType::Scalar Scalar;
303 typedef JacobiSVD<MatrixType, Options> SVDType;
306 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
307 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
308 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
309 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
310 MatrixOptions = internal::traits<MatrixType>::Options,
311 WorkspaceSize = internal::traits<SVDType>::MatrixVColsAtCompileTime,
312 MaxWorkspaceSize = internal::traits<SVDType>::MatrixVMaxColsAtCompileTime
315 typedef Matrix<Scalar, WorkspaceSize, 1, ColMajor, MaxWorkspaceSize, 1> WorkspaceType;
317 typedef typename internal::make_proper_matrix_type<Scalar, ColsAtCompileTime, RowsAtCompileTime, MatrixOptions,
318 MaxColsAtCompileTime, MaxRowsAtCompileTime>::type
319 TransposeTypeWithSameStorageOrder;
321 void allocate(
const SVDType& svd) {
322 if (svd.cols() != m_qr.rows() || svd.rows() != m_qr.cols()) {
323 internal::destroy_at(&m_qr);
324 internal::construct_at(&m_qr, svd.cols(), svd.rows());
326 if (svd.m_computeFullV)
327 m_workspace.resize(svd.cols());
328 else if (svd.m_computeThinV)
329 m_workspace.resize(svd.rows());
332 template <
typename Xpr>
333 bool run(SVDType& svd,
const Xpr& matrix) {
334 if (matrix.cols() > matrix.rows()) {
335 m_qr.compute(matrix.adjoint());
338 m_qr.matrixQR().block(0, 0, matrix.rows(), matrix.rows()).template triangularView<Upper>().adjoint();
339 if (svd.m_computeFullV)
340 m_qr.householderQ().evalTo(svd.m_matrixV, m_workspace);
341 else if (svd.m_computeThinV) {
342 svd.m_matrixV.setIdentity(matrix.cols(), matrix.rows());
343 m_qr.householderQ().applyThisOnTheLeft(svd.m_matrixV, m_workspace);
345 if (svd.computeU()) svd.m_matrixU.setIdentity(matrix.rows(), matrix.rows());
352 typedef HouseholderQR<TransposeTypeWithSameStorageOrder> QRType;
354 WorkspaceType m_workspace;
362template <
typename MatrixType,
int Options>
363struct svd_precondition_2x2_block_to_be_real<MatrixType, Options, false> {
364 typedef JacobiSVD<MatrixType, Options> SVD;
365 typedef typename MatrixType::RealScalar RealScalar;
366 static bool run(
typename SVD::WorkMatrixType&, SVD&,
Index,
Index, RealScalar&) {
return true; }
369template <
typename MatrixType,
int Options>
370struct svd_precondition_2x2_block_to_be_real<MatrixType, Options, true> {
371 typedef JacobiSVD<MatrixType, Options> SVD;
372 typedef typename MatrixType::Scalar Scalar;
373 typedef typename MatrixType::RealScalar RealScalar;
374 static bool run(
typename SVD::WorkMatrixType& work_matrix, SVD& svd,
Index p,
Index q, RealScalar& maxDiagEntry) {
378 JacobiRotation<Scalar> rot;
379 RealScalar n =
sqrt(numext::abs2(work_matrix.coeff(p, p)) + numext::abs2(work_matrix.coeff(q, p)));
381 const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
382 const RealScalar precision = NumTraits<Scalar>::epsilon();
384 if (numext::is_exactly_zero(n)) {
386 work_matrix.coeffRef(p, p) = work_matrix.coeffRef(q, p) = Scalar(0);
388 if (
abs(numext::imag(work_matrix.coeff(p, q))) > considerAsZero) {
391 z =
abs(work_matrix.coeff(p, q)) / work_matrix.coeff(p, q);
392 work_matrix.row(p) *= z;
393 if (svd.computeU()) svd.m_matrixU.col(p) *=
conj(z);
395 if (
abs(numext::imag(work_matrix.coeff(q, q))) > considerAsZero) {
396 z =
abs(work_matrix.coeff(q, q)) / work_matrix.coeff(q, q);
397 work_matrix.row(q) *= z;
398 if (svd.computeU()) svd.m_matrixU.col(q) *=
conj(z);
402 rot.c() =
conj(work_matrix.coeff(p, p)) / n;
403 rot.s() = work_matrix.coeff(q, p) / n;
404 work_matrix.applyOnTheLeft(p, q, rot);
405 if (svd.computeU()) svd.m_matrixU.applyOnTheRight(p, q, rot.adjoint());
406 if (
abs(numext::imag(work_matrix.coeff(p, q))) > considerAsZero) {
407 z =
abs(work_matrix.coeff(p, q)) / work_matrix.coeff(p, q);
408 work_matrix.col(q) *= z;
409 if (svd.computeV()) svd.m_matrixV.col(q) *= z;
411 if (
abs(numext::imag(work_matrix.coeff(q, q))) > considerAsZero) {
412 z =
abs(work_matrix.coeff(q, q)) / work_matrix.coeff(q, q);
413 work_matrix.row(q) *= z;
414 if (svd.computeU()) svd.m_matrixU.col(q) *=
conj(z);
419 maxDiagEntry = numext::maxi<RealScalar>(
420 maxDiagEntry, numext::maxi<RealScalar>(
abs(work_matrix.coeff(p, p)),
abs(work_matrix.coeff(q, q))));
422 RealScalar threshold = numext::maxi<RealScalar>(considerAsZero, precision * maxDiagEntry);
423 return abs(work_matrix.coeff(p, q)) > threshold ||
abs(work_matrix.coeff(q, p)) > threshold;
427template <
typename MatrixType_,
int Options>
428struct traits<JacobiSVD<MatrixType_, Options> > : svd_traits<MatrixType_, Options> {
429 typedef MatrixType_ MatrixType;
499template <
typename MatrixType_,
int Options_>
504 typedef MatrixType_ MatrixType;
505 typedef typename Base::Scalar Scalar;
506 typedef typename Base::RealScalar RealScalar;
509 QRPreconditioner = internal::get_qr_preconditioner(Options),
510 RowsAtCompileTime = Base::RowsAtCompileTime,
511 ColsAtCompileTime = Base::ColsAtCompileTime,
512 DiagSizeAtCompileTime = Base::DiagSizeAtCompileTime,
513 MaxRowsAtCompileTime = Base::MaxRowsAtCompileTime,
514 MaxColsAtCompileTime = Base::MaxColsAtCompileTime,
515 MaxDiagSizeAtCompileTime = Base::MaxDiagSizeAtCompileTime,
516 MatrixOptions = Base::MatrixOptions
519 typedef typename Base::MatrixUType MatrixUType;
520 typedef typename Base::MatrixVType MatrixVType;
521 typedef typename Base::SingularValuesType SingularValuesType;
522 typedef Matrix<Scalar, DiagSizeAtCompileTime, DiagSizeAtCompileTime, MatrixOptions, MaxDiagSizeAtCompileTime,
523 MaxDiagSizeAtCompileTime>
560 internal::check_svd_options_assertions<MatrixType, Options>(computationOptions, rows, cols);
561 allocate(rows, cols, computationOptions);
569 template <
typename Derived>
571 compute_impl(matrix, internal::get_computation_options(Options));
587 template <
typename Derived>
589 internal::check_svd_options_assertions<MatrixBase<Derived>, Options>(computationOptions, matrix.rows(),
591 compute_impl(matrix, computationOptions);
599 template <
typename Derived>
601 return compute_impl(matrix, m_computationOptions);
613 template <
typename Derived>
616 internal::check_svd_options_assertions<MatrixBase<Derived>, Options>(m_computationOptions, matrix.rows(),
618 return compute_impl(matrix, computationOptions);
624 using Base::diagSize;
628 void allocate(
Index rows_,
Index cols_,
unsigned int computationOptions) {
629 if (Base::allocate(rows_, cols_, computationOptions))
return;
632 "JacobiSVD: can't compute thin U or thin V with the FullPivHouseholderQR preconditioner. "
633 "Use the ColPivHouseholderQR preconditioner instead.");
635 m_workMatrix.
resize(diagSize(), diagSize());
636 if (
cols() >
rows()) m_qr_precond_morecols.allocate(*
this);
637 if (
rows() >
cols()) m_qr_precond_morerows.allocate(*
this);
641 template <
typename Derived>
642 JacobiSVD& compute_impl(
const MatrixBase<Derived>& matrix,
unsigned int computationOptions);
645 using Base::m_computationOptions;
646 using Base::m_computeFullU;
647 using Base::m_computeFullV;
648 using Base::m_computeThinU;
649 using Base::m_computeThinV;
651 using Base::m_isAllocated;
652 using Base::m_isInitialized;
653 using Base::m_matrixU;
654 using Base::m_matrixV;
655 using Base::m_nonzeroSingularValues;
656 using Base::m_prescribedThreshold;
657 using Base::m_singularValues;
658 using Base::m_usePrescribedThreshold;
659 using Base::ShouldComputeThinU;
660 using Base::ShouldComputeThinV;
664 "JacobiSVD: can't compute thin U or thin V with the FullPivHouseholderQR preconditioner. "
665 "Use the ColPivHouseholderQR preconditioner instead.")
667 template <typename MatrixType__,
int Options__,
bool IsComplex_>
668 friend struct internal::svd_precondition_2x2_block_to_be_real;
669 template <typename MatrixType__,
int Options__,
int QRPreconditioner_,
int Case_,
bool DoAnything_>
670 friend struct internal::qr_preconditioner_impl;
672 internal::qr_preconditioner_impl<MatrixType, Options, QRPreconditioner, internal::PreconditionIfMoreColsThanRows>
673 m_qr_precond_morecols;
674 internal::qr_preconditioner_impl<MatrixType, Options, QRPreconditioner, internal::PreconditionIfMoreRowsThanCols>
675 m_qr_precond_morerows;
676 WorkMatrixType m_workMatrix;
679template <typename MatrixType,
int Options>
680template <typename Derived>
681JacobiSVD<MatrixType, Options>&
JacobiSVD<MatrixType, Options>::compute_impl(const MatrixBase<Derived>& matrix,
682 unsigned int computationOptions) {
683 EIGEN_STATIC_ASSERT_SAME_MATRIX_SIZE(Derived, MatrixType);
684 EIGEN_STATIC_ASSERT((std::is_same<typename Derived::Scalar, typename MatrixType::Scalar>::value),
685 Input matrix must have the same Scalar type as the BDCSVD
object.);
689 allocate(matrix.rows(), matrix.cols(), computationOptions);
693 const RealScalar precision = RealScalar(2) * NumTraits<Scalar>::epsilon();
696 const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
699 RealScalar scale = matrix.cwiseAbs().template maxCoeff<PropagateNaN>();
700 if (!(numext::isfinite)(scale)) {
701 m_isInitialized =
true;
703 m_nonzeroSingularValues = 0;
706 if (numext::is_exactly_zero(scale)) scale = RealScalar(1);
711 m_qr_precond_morecols.run(*
this, matrix / scale);
712 m_qr_precond_morerows.run(*
this, matrix / scale);
715 matrix.template topLeftCorner<DiagSizeAtCompileTime, DiagSizeAtCompileTime>(diagSize(), diagSize()) / scale;
716 if (m_computeFullU) m_matrixU.setIdentity(
rows(),
rows());
717 if (m_computeThinU) m_matrixU.setIdentity(
rows(), diagSize());
718 if (m_computeFullV) m_matrixV.setIdentity(
cols(),
cols());
719 if (m_computeThinV) m_matrixV.setIdentity(
cols(), diagSize());
723 RealScalar maxDiagEntry = m_workMatrix.cwiseAbs().diagonal().maxCoeff();
725 bool finished =
false;
731 for (
Index p = 1; p < diagSize(); ++p) {
732 for (
Index q = 0; q < p; ++q) {
736 RealScalar
threshold = numext::maxi<RealScalar>(considerAsZero, precision * maxDiagEntry);
741 if (internal::svd_precondition_2x2_block_to_be_real<MatrixType, Options>::run(m_workMatrix, *
this, p, q,
743 JacobiRotation<RealScalar> j_left, j_right;
744 internal::real_2x2_jacobi_svd(m_workMatrix, p, q, &j_left, &j_right);
747 m_workMatrix.applyOnTheLeft(p, q, j_left);
748 if (
computeU()) m_matrixU.applyOnTheRight(p, q, j_left.transpose());
750 m_workMatrix.applyOnTheRight(p, q, j_right);
751 if (
computeV()) m_matrixV.applyOnTheRight(p, q, j_right);
754 maxDiagEntry = numext::maxi<RealScalar>(
755 maxDiagEntry, numext::maxi<RealScalar>(
abs(m_workMatrix.coeff(p, p)),
abs(m_workMatrix.coeff(q, q))));
765 for (
Index i = 0; i < diagSize(); ++i) {
769 if (NumTraits<Scalar>::IsComplex &&
abs(numext::imag(m_workMatrix.coeff(i, i))) > considerAsZero) {
770 RealScalar a =
abs(m_workMatrix.coeff(i, i));
771 m_singularValues.coeffRef(i) =
abs(a);
772 if (
computeU()) m_matrixU.col(i) *= m_workMatrix.coeff(i, i) / a;
775 RealScalar a = numext::real(m_workMatrix.coeff(i, i));
776 m_singularValues.coeffRef(i) =
abs(a);
777 if (
computeU() && (a < RealScalar(0))) m_matrixU.col(i) = -m_matrixU.col(i);
781 m_singularValues *= scale;
785 m_nonzeroSingularValues = diagSize();
786 for (
Index i = 0; i < diagSize(); i++) {
788 RealScalar maxRemainingSingularValue = m_singularValues.tail(diagSize() - i).maxCoeff(&pos);
789 if (numext::is_exactly_zero(maxRemainingSingularValue)) {
790 m_nonzeroSingularValues = i;
795 std::swap(m_singularValues.coeffRef(i), m_singularValues.coeffRef(pos));
796 if (
computeU()) m_matrixU.col(pos).swap(m_matrixU.col(i));
797 if (
computeV()) m_matrixV.col(pos).swap(m_matrixV.col(i));
801 m_isInitialized =
true;
812template <
typename Derived>
813template <
int Options>
818template <
typename Derived>
819template <
int Options>
821 unsigned int computationOptions)
const {
Two-sided Jacobi SVD decomposition of a rectangular matrix.
Definition JacobiSVD.h:500
JacobiSVD()
Default Constructor.
Definition JacobiSVD.h:531
bool computeV() const
Definition SVDBase.h:275
JacobiSVD & compute(const MatrixBase< Derived > &matrix)
Method performing the decomposition of given matrix. Computes Thin/Full unitaries U/V if specified us...
Definition JacobiSVD.h:600
bool computeU() const
Definition SVDBase.h:273
JacobiSVD(Index rows, Index cols)
Default Constructor with memory preallocation.
Definition JacobiSVD.h:540
JacobiSVD(const MatrixBase< Derived > &matrix)
Constructor performing the decomposition of given matrix, using the custom options specified with the...
Definition JacobiSVD.h:570
JacobiSVD(const MatrixBase< Derived > &matrix, unsigned int computationOptions)
Constructor performing the decomposition of given matrix using specified options for computing unitar...
Definition JacobiSVD.h:588
Base class for all dense matrices, vectors, and expressions.
Definition MatrixBase.h:52
JacobiSVD< typename MatrixBase< Derived >::PlainObject, Options > jacobiSvd() const
Definition JacobiSVD.h:814
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:186
constexpr void resize(Index rows, Index cols)
Definition PlainObjectBase.h:268
Index rank() const
Definition SVDBase.h:217
bool computeV() const
Definition SVDBase.h:275
bool computeU() const
Definition SVDBase.h:273
RealScalar threshold() const
Definition SVDBase.h:265
SVDBase()
Definition SVDBase.h:349
@ NoQRPreconditioner
Definition Constants.h:423
@ HouseholderQRPreconditioner
Definition Constants.h:425
@ ColPivHouseholderQRPreconditioner
Definition Constants.h:421
@ FullPivHouseholderQRPreconditioner
Definition Constants.h:427
@ InvalidInput
Definition Constants.h:447
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
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_conjugate_op< typename Derived::Scalar >, const Derived > conj(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:82
const int Dynamic
Definition Constants.h:25
constexpr Index cols() const noexcept
Definition EigenBase.h:61
Eigen::Index Index
Definition EigenBase.h:43
constexpr Index rows() const noexcept
Definition EigenBase.h:59