11#ifndef EIGEN_JACOBISVD_H
12#define EIGEN_JACOBISVD_H
19template<
typename MatrixType,
int QRPreconditioner,
20 bool IsComplex = NumTraits<typename MatrixType::Scalar>::IsComplex>
30enum { PreconditionIfMoreColsThanRows, PreconditionIfMoreRowsThanCols };
32template<
typename MatrixType,
int QRPreconditioner,
int Case>
35 enum { a = MatrixType::RowsAtCompileTime !=
Dynamic &&
36 MatrixType::ColsAtCompileTime !=
Dynamic &&
37 MatrixType::ColsAtCompileTime <= MatrixType::RowsAtCompileTime,
38 b = MatrixType::RowsAtCompileTime !=
Dynamic &&
39 MatrixType::ColsAtCompileTime !=
Dynamic &&
40 MatrixType::RowsAtCompileTime <= MatrixType::ColsAtCompileTime,
42 (
Case == PreconditionIfMoreColsThanRows &&
bool(a)) ||
43 (
Case == PreconditionIfMoreRowsThanCols &&
bool(b)) )
51template<
typename MatrixType,
int QRPreconditioner,
int Case>
64template<
typename MatrixType>
68 typedef typename MatrixType::Scalar Scalar;
71 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
72 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime
78 if (
svd.rows() != m_qr.rows() ||
svd.cols() != m_qr.cols())
83 if (
svd.m_computeFullU) m_workspace.resize(
svd.rows());
88 if(matrix.rows() > matrix.cols())
92 if(
svd.m_computeFullU) m_qr.matrixQ().evalTo(
svd.m_matrixU, m_workspace);
93 if(
svd.computeV())
svd.m_matrixV = m_qr.colsPermutation();
104template<
typename MatrixType>
108 typedef typename MatrixType::Scalar Scalar;
111 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
112 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
113 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
114 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
115 Options = MatrixType::Options
119 Scalar, ColsAtCompileTime, RowsAtCompileTime, Options, MaxColsAtCompileTime, MaxRowsAtCompileTime
124 if (
svd.cols() != m_qr.rows() ||
svd.rows() != m_qr.cols())
129 m_adjoint.resize(
svd.cols(),
svd.rows());
130 if (
svd.m_computeFullV) m_workspace.resize(
svd.cols());
135 if(matrix.cols() > matrix.rows())
137 m_adjoint = matrix.adjoint();
138 m_qr.compute(m_adjoint);
140 if(
svd.m_computeFullV) m_qr.matrixQ().evalTo(
svd.m_matrixV, m_workspace);
141 if(
svd.computeU())
svd.m_matrixU = m_qr.colsPermutation();
155template<
typename MatrixType>
161 if (
svd.rows() != m_qr.rows() ||
svd.cols() != m_qr.cols())
166 if (
svd.m_computeFullU) m_workspace.resize(
svd.rows());
167 else if (
svd.m_computeThinU) m_workspace.resize(
svd.cols());
172 if(matrix.rows() > matrix.cols())
174 m_qr.compute(matrix);
176 if(
svd.m_computeFullU) m_qr.householderQ().evalTo(
svd.m_matrixU, m_workspace);
177 else if(
svd.m_computeThinU)
180 m_qr.householderQ().applyThisOnTheLeft(
svd.m_matrixU, m_workspace);
182 if(
svd.computeV())
svd.m_matrixV = m_qr.colsPermutation();
194template<
typename MatrixType>
198 typedef typename MatrixType::Scalar Scalar;
201 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
202 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
203 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
204 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
205 Options = MatrixType::Options
209 Scalar, ColsAtCompileTime, RowsAtCompileTime, Options, MaxColsAtCompileTime, MaxRowsAtCompileTime
214 if (
svd.cols() != m_qr.rows() ||
svd.rows() != m_qr.cols())
219 if (
svd.m_computeFullV) m_workspace.resize(
svd.cols());
220 else if (
svd.m_computeThinV) m_workspace.resize(
svd.rows());
221 m_adjoint.resize(
svd.cols(),
svd.rows());
226 if(matrix.cols() > matrix.rows())
228 m_adjoint = matrix.adjoint();
229 m_qr.compute(m_adjoint);
232 if(
svd.m_computeFullV) m_qr.householderQ().evalTo(
svd.m_matrixV, m_workspace);
233 else if(
svd.m_computeThinV)
236 m_qr.householderQ().applyThisOnTheLeft(
svd.m_matrixV, m_workspace);
238 if(
svd.computeU())
svd.m_matrixU = m_qr.colsPermutation();
253template<
typename MatrixType>
259 if (
svd.rows() != m_qr.rows() ||
svd.cols() != m_qr.cols())
264 if (
svd.m_computeFullU) m_workspace.resize(
svd.rows());
265 else if (
svd.m_computeThinU) m_workspace.resize(
svd.cols());
270 if(matrix.rows() > matrix.cols())
272 m_qr.compute(matrix);
274 if(
svd.m_computeFullU) m_qr.householderQ().evalTo(
svd.m_matrixU, m_workspace);
275 else if(
svd.m_computeThinU)
278 m_qr.householderQ().applyThisOnTheLeft(
svd.m_matrixU, m_workspace);
291template<
typename MatrixType>
295 typedef typename MatrixType::Scalar Scalar;
298 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
299 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
300 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
301 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
302 Options = MatrixType::Options
306 Scalar, ColsAtCompileTime, RowsAtCompileTime, Options, MaxColsAtCompileTime, MaxRowsAtCompileTime
311 if (
svd.cols() != m_qr.rows() ||
svd.rows() != m_qr.cols())
316 if (
svd.m_computeFullV) m_workspace.resize(
svd.cols());
317 else if (
svd.m_computeThinV) m_workspace.resize(
svd.rows());
318 m_adjoint.resize(
svd.cols(),
svd.rows());
323 if(matrix.cols() > matrix.rows())
325 m_adjoint = matrix.adjoint();
326 m_qr.compute(m_adjoint);
329 if(
svd.m_computeFullV) m_qr.householderQ().evalTo(
svd.m_matrixV, m_workspace);
330 else if(
svd.m_computeThinV)
333 m_qr.householderQ().applyThisOnTheLeft(
svd.m_matrixV, m_workspace);
353template<
typename MatrixType,
int QRPreconditioner>
357 typedef typename MatrixType::RealScalar RealScalar;
358 static bool run(
typename SVD::WorkMatrixType&,
SVD&,
Index,
Index, RealScalar&) {
return true; }
361template<
typename MatrixType,
int QRPreconditioner>
365 typedef typename MatrixType::Scalar Scalar;
366 typedef typename MatrixType::RealScalar
RealScalar;
388 if(
svd.computeU())
svd.m_matrixU.col(p) *=
conj(z);
394 if(
svd.computeU())
svd.m_matrixU.col(
q) *=
conj(z);
408 if(
svd.computeV())
svd.m_matrixV.col(
q) *= z;
414 if(
svd.computeU())
svd.m_matrixU.col(
q) *=
conj(z);
426template<
typename _MatrixType,
int QRPreconditioner>
430 typedef _MatrixType MatrixType;
488template<
typename _MatrixType,
int QRPreconditioner>
class JacobiSVD
489 :
public SVDBase<JacobiSVD<_MatrixType,QRPreconditioner> >
494 typedef _MatrixType MatrixType;
495 typedef typename MatrixType::Scalar Scalar;
498 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
499 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
500 DiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_DYNAMIC(RowsAtCompileTime,ColsAtCompileTime),
501 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
502 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
503 MaxDiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(MaxRowsAtCompileTime,MaxColsAtCompileTime),
504 MatrixOptions = MatrixType::Options
513 typedef Matrix<Scalar, DiagSizeAtCompileTime, DiagSizeAtCompileTime,
514 MatrixOptions, MaxDiagSizeAtCompileTime, MaxDiagSizeAtCompileTime>
572 return compute(matrix, m_computationOptions);
585 using Base::m_matrixU;
586 using Base::m_matrixV;
587 using Base::m_singularValues;
589 using Base::m_isInitialized;
590 using Base::m_isAllocated;
591 using Base::m_usePrescribedThreshold;
592 using Base::m_computeFullU;
593 using Base::m_computeThinU;
594 using Base::m_computeFullV;
595 using Base::m_computeThinV;
596 using Base::m_computationOptions;
597 using Base::m_nonzeroSingularValues;
600 using Base::m_diagSize;
601 using Base::m_prescribedThreshold;
602 WorkMatrixType m_workMatrix;
604 template<
typename __MatrixType,
int _QRPreconditioner,
bool _IsComplex>
606 template<
typename __MatrixType,
int _QRPreconditioner,
int _Case,
bool _DoAnything>
611 MatrixType m_scaledMatrix;
614template<
typename MatrixType,
int QRPreconditioner>
615void JacobiSVD<MatrixType, QRPreconditioner>::allocate(
Eigen::Index rows,
Eigen::Index cols,
unsigned int computationOptions)
617 eigen_assert(rows >= 0 && cols >= 0);
622 computationOptions == m_computationOptions)
630 m_isInitialized =
false;
631 m_isAllocated =
true;
632 m_computationOptions = computationOptions;
633 m_computeFullU = (computationOptions &
ComputeFullU) != 0;
634 m_computeThinU = (computationOptions &
ComputeThinU) != 0;
635 m_computeFullV = (computationOptions &
ComputeFullV) != 0;
636 m_computeThinV = (computationOptions &
ComputeThinV) != 0;
637 eigen_assert(!(m_computeFullU && m_computeThinU) &&
"JacobiSVD: you can't ask for both full and thin U");
638 eigen_assert(!(m_computeFullV && m_computeThinV) &&
"JacobiSVD: you can't ask for both full and thin V");
639 eigen_assert(EIGEN_IMPLIES(m_computeThinU || m_computeThinV, MatrixType::ColsAtCompileTime==
Dynamic) &&
640 "JacobiSVD: thin U and V are only available when your matrix has a dynamic number of columns.");
643 eigen_assert(!(m_computeThinU || m_computeThinV) &&
644 "JacobiSVD: can't compute thin U or thin V with the FullPivHouseholderQR preconditioner. "
645 "Use the ColPivHouseholderQR preconditioner instead.");
647 m_diagSize = (std::min)(m_rows, m_cols);
648 m_singularValues.resize(m_diagSize);
650 m_matrixU.resize(m_rows, m_computeFullU ? m_rows
651 : m_computeThinU ? m_diagSize
654 m_matrixV.resize(m_cols, m_computeFullV ? m_cols
655 : m_computeThinV ? m_diagSize
657 m_workMatrix.resize(m_diagSize, m_diagSize);
659 if(m_cols>m_rows) m_qr_precond_morecols.allocate(*
this);
660 if(m_rows>m_cols) m_qr_precond_morerows.allocate(*
this);
661 if(m_rows!=m_cols) m_scaledMatrix.resize(rows,cols);
664template<
typename MatrixType,
int QRPreconditioner>
665JacobiSVD<MatrixType, QRPreconditioner>&
680 if (!(numext::isfinite)(scale)) {
681 m_isInitialized =
true;
683 m_nonzeroSingularValues = 0;
692 m_scaledMatrix = matrix / scale;
693 m_qr_precond_morecols.run(*
this, m_scaledMatrix);
694 m_qr_precond_morerows.run(*
this, m_scaledMatrix);
698 m_workMatrix = matrix.block(0,0,m_diagSize,m_diagSize) / scale;
699 if(m_computeFullU) m_matrixU.
setIdentity(m_rows,m_rows);
700 if(m_computeThinU) m_matrixU.setIdentity(m_rows,m_diagSize);
701 if(m_computeFullV) m_matrixV.setIdentity(m_cols,m_cols);
702 if(m_computeThinV) m_matrixV.setIdentity(m_cols, m_diagSize);
708 bool finished =
false;
715 for(
Index p = 1; p < m_diagSize; ++p)
723 if(abs(m_workMatrix.coeff(p,
q))>threshold || abs(m_workMatrix.coeff(
q,p)) > threshold)
731 internal::real_2x2_jacobi_svd(m_workMatrix, p,
q, &
j_left, &
j_right);
734 m_workMatrix.applyOnTheLeft(p,
q,
j_left);
737 m_workMatrix.applyOnTheRight(p,
q,
j_right);
738 if(computeV()) m_matrixV.applyOnTheRight(p,
q,
j_right);
741 maxDiagEntry = numext::maxi<RealScalar>(
maxDiagEntry,numext::maxi<RealScalar>(abs(m_workMatrix.coeff(p,p)), abs(m_workMatrix.coeff(
q,
q))));
750 for(
Index i = 0; i < m_diagSize; ++i)
758 m_singularValues.coeffRef(i) = abs(a);
759 if(computeU()) m_matrixU.col(i) *= m_workMatrix.coeff(i,i)/a;
764 RealScalar a = numext::real(m_workMatrix.coeff(i,i));
765 m_singularValues.coeffRef(i) = abs(a);
766 if(computeU() && (a<
RealScalar(0))) m_matrixU.col(i) = -m_matrixU.col(i);
770 m_singularValues *= scale;
774 m_nonzeroSingularValues = m_diagSize;
775 for(
Index i = 0; i < m_diagSize; i++)
781 m_nonzeroSingularValues = i;
787 std::swap(m_singularValues.coeffRef(i), m_singularValues.coeffRef(pos));
788 if(computeU()) m_matrixU.col(pos).swap(m_matrixU.col(i));
789 if(computeV()) m_matrixV.col(pos).swap(m_matrixV.col(i));
793 m_isInitialized =
true;
804template<
typename Derived>
EIGEN_DEVICE_FUNC TransposeReturnType transpose()
Definition Transpose.h:182
Two-sided Jacobi SVD decomposition of a rectangular matrix.
Definition JacobiSVD.h:490
JacobiSVD & compute(const MatrixType &matrix)
Method performing the decomposition of given matrix using current options.
Definition JacobiSVD.h:570
JacobiSVD()
Default Constructor.
Definition JacobiSVD.h:522
JacobiSVD(Index rows, Index cols, unsigned int computationOptions=0)
Default Constructor with memory preallocation.
Definition JacobiSVD.h:532
JacobiSVD & compute(const MatrixType &matrix, unsigned int computationOptions)
Method performing the decomposition of given matrix using custom options.
Definition JacobiSVD.h:666
JacobiSVD(const MatrixType &matrix, unsigned int computationOptions=0)
Constructor performing the decomposition of given matrix.
Definition JacobiSVD.h:547
Base class for all dense matrices, vectors, and expressions.
Definition MatrixBase.h:50
void applyOnTheLeft(const EigenBase< OtherDerived > &other)
replaces *this by other * *this.
Definition MatrixBase.h:534
void applyOnTheRight(const EigenBase< OtherDerived > &other)
replaces *this by *this * other.
Definition MatrixBase.h:522
EIGEN_DEVICE_FUNC const AdjointReturnType adjoint() const
Definition Transpose.h:221
@ ColsAtCompileTime
The number of columns at compile-time.
Definition DenseBase.h:106
@ RowsAtCompileTime
The number of rows at compile-time.
Definition DenseBase.h:100
EIGEN_DEVICE_FUNC Derived & setIdentity()
Writes the identity expression (not necessarily square) into *this.
Definition CwiseNullaryOp.h:873
EIGEN_DEVICE_FUNC DiagonalReturnType diagonal()
Definition Diagonal.h:187
Base class of SVD algorithms.
Definition SVDBase.h:64
Index rank() const
Definition SVDBase.h:148
bool computeV() const
Definition SVDBase.h:210
Eigen::Index Index
Definition SVDBase.h:74
bool computeU() const
Definition SVDBase.h:208
Definition XprHelper.h:258
@ NoQRPreconditioner
Do not specify what is to be done if the SVD of a non-square matrix is asked for.
Definition Constants.h:425
@ HouseholderQRPreconditioner
Use a QR decomposition without pivoting as the first step.
Definition Constants.h:427
@ ColPivHouseholderQRPreconditioner
Use a QR decomposition with column pivoting as the first step.
Definition Constants.h:429
@ FullPivHouseholderQRPreconditioner
Use a QR decomposition with full pivoting as the first step.
Definition Constants.h:431
@ InvalidInput
The inputs are invalid, or the algorithm has been improperly called.
Definition Constants.h:449
@ Success
Computation was successful.
Definition Constants.h:442
@ ComputeFullV
Used in JacobiSVD to indicate that the square matrix V is to be computed.
Definition Constants.h:397
@ ComputeThinV
Used in JacobiSVD to indicate that the thin matrix V is to be computed.
Definition Constants.h:399
@ ComputeFullU
Used in JacobiSVD to indicate that the square matrix U is to be computed.
Definition Constants.h:393
@ ComputeThinU
Used in JacobiSVD to indicate that the thin matrix U is to be computed.
Definition Constants.h:395
Namespace containing all symbols from the Eigen library.
Definition LDLT.h:16
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:74
const int Dynamic
This value means that a positive quantity (e.g., a size) is not known at compile-time,...
Definition Constants.h:22
NLOHMANN_BASIC_JSON_TPL_DECLARATION void swap(nlohmann::NLOHMANN_BASIC_JSON_TPL &j1, nlohmann::NLOHMANN_BASIC_JSON_TPL &j2) noexcept(//NOLINT(readability-inconsistent-declaration-parameter-name, cert-dcl58-cpp) is_nothrow_move_constructible< nlohmann::NLOHMANN_BASIC_JSON_TPL >::value &&//NOLINT(misc-redundant-expression) is_nothrow_move_assignable< nlohmann::NLOHMANN_BASIC_JSON_TPL >::value)
exchanges the values of two JSON objects
Definition json.hpp:24418
Holds information about the various numeric (i.e.
Definition NumTraits.h:236
Definition JacobiSVD.h:49
Definition JacobiSVD.h:34
Definition JacobiSVD.h:21
Definition ForwardDeclarations.h:17