|
enum | : int {
Options = Options_
, QRPreconditioner = internal::get_qr_preconditioner(Options)
, RowsAtCompileTime = Base::RowsAtCompileTime
, ColsAtCompileTime = Base::ColsAtCompileTime
,
DiagSizeAtCompileTime = Base::DiagSizeAtCompileTime
, MaxRowsAtCompileTime = Base::MaxRowsAtCompileTime
, MaxColsAtCompileTime = Base::MaxColsAtCompileTime
, MaxDiagSizeAtCompileTime = Base::MaxDiagSizeAtCompileTime
,
MatrixOptions = Base::MatrixOptions
} |
|
typedef MatrixType_ | MatrixType |
|
typedef Base::Scalar | Scalar |
|
typedef Base::RealScalar | RealScalar |
|
typedef Base::MatrixUType | MatrixUType |
|
typedef Base::MatrixVType | MatrixVType |
|
typedef Base::SingularValuesType | SingularValuesType |
|
typedef Matrix< Scalar, DiagSizeAtCompileTime, DiagSizeAtCompileTime, MatrixOptions, MaxDiagSizeAtCompileTime, MaxDiagSizeAtCompileTime > | WorkMatrixType |
|
enum | |
|
typedef internal::traits< JacobiSVD< MatrixType_, Options_ > >::MatrixType | MatrixType |
|
typedef MatrixType::Scalar | Scalar |
|
typedef NumTraits< typename MatrixType::Scalar >::Real | RealScalar |
|
typedef Eigen::internal::traits< SVDBase >::StorageIndex | StorageIndex |
|
typedef internal::make_proper_matrix_type< Scalar, RowsAtCompileTime, MatrixUColsAtCompileTime, MatrixOptions, MaxRowsAtCompileTime, MatrixUMaxColsAtCompileTime >::type | MatrixUType |
|
typedef internal::make_proper_matrix_type< Scalar, ColsAtCompileTime, MatrixVColsAtCompileTime, MatrixOptions, MaxColsAtCompileTime, MatrixVMaxColsAtCompileTime >::type | MatrixVType |
|
typedef internal::plain_diag_type< MatrixType, RealScalar >::type | SingularValuesType |
|
enum | {
RowsAtCompileTime = internal::traits<Derived>::RowsAtCompileTime
, ColsAtCompileTime = internal::traits<Derived>::ColsAtCompileTime
, SizeAtCompileTime = (internal::size_of_xpr_at_compile_time<Derived>::ret)
, MaxRowsAtCompileTime = internal::traits<Derived>::MaxRowsAtCompileTime
,
MaxColsAtCompileTime = internal::traits<Derived>::MaxColsAtCompileTime
, MaxSizeAtCompileTime
, IsVectorAtCompileTime
, NumDimensions
} |
|
typedef EigenBase< Derived > | Base |
|
typedef internal::traits< Derived >::Scalar | Scalar |
|
typedef Scalar | CoeffReturnType |
|
typedef Transpose< const Derived > | ConstTransposeReturnType |
|
typedef std::conditional_t< NumTraits< Scalar >::IsComplex, CwiseUnaryOp< internal::scalar_conjugate_op< Scalar >, const ConstTransposeReturnType >, const ConstTransposeReturnType > | AdjointReturnType |
|
typedef Eigen::Index | Index |
| The interface type of indices. More...
|
|
typedef internal::traits< Derived >::StorageKind | StorageKind |
|
|
| JacobiSVD () |
| Default Constructor. More...
|
|
| JacobiSVD (Index rows, Index cols) |
| Default Constructor with memory preallocation. More...
|
|
EIGEN_DEPRECATED | JacobiSVD (Index rows, Index cols, unsigned int computationOptions) |
| Default Constructor with memory preallocation. More...
|
|
| JacobiSVD (const MatrixType &matrix) |
| Constructor performing the decomposition of given matrix, using the custom options specified with the Options template parameter. More...
|
|
| JacobiSVD (const MatrixType &matrix, unsigned int computationOptions) |
| Constructor performing the decomposition of given matrix using specified options for computing unitaries. More...
|
|
JacobiSVD & | compute (const MatrixType &matrix) |
| Method performing the decomposition of given matrix. Computes Thin/Full unitaries U/V if specified using the Options template parameter or the class constructor. More...
|
|
EIGEN_DEPRECATED JacobiSVD & | compute (const MatrixType &matrix, unsigned int computationOptions) |
| Method performing the decomposition of given matrix, as specified by the computationOptions parameter. More...
|
|
void | allocate (Index rows_, Index cols_, unsigned int computationOptions) |
|
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index | cols () const EIGEN_NOEXCEPT |
|
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index | rows () const EIGEN_NOEXCEPT |
|
JacobiSVD< MatrixType_, Options_ > & | derived () |
|
const JacobiSVD< MatrixType_, Options_ > & | derived () const |
|
const MatrixUType & | matrixU () const |
|
const MatrixVType & | matrixV () const |
|
const SingularValuesType & | singularValues () const |
|
Index | nonzeroSingularValues () const |
|
Index | rank () const |
|
JacobiSVD< MatrixType_, Options_ > & | setThreshold (const RealScalar &threshold) |
|
JacobiSVD< MatrixType_, Options_ > & | setThreshold (Default_t) |
|
RealScalar | threshold () const |
|
bool | computeU () const |
|
bool | computeV () const |
|
Index | rows () const |
|
Index | cols () const |
|
Index | diagSize () const |
|
EIGEN_DEVICE_FUNC ComputationInfo | info () const |
| Reports whether previous computation was successful. More...
|
|
void | _solve_impl (const RhsType &rhs, DstType &dst) const |
|
void | _solve_impl_transposed (const RhsType &rhs, DstType &dst) const |
|
| SolverBase () |
|
| ~SolverBase () |
|
template<typename Rhs > |
const Solve< Derived, Rhs > | solve (const MatrixBase< Rhs > &b) const |
|
const ConstTransposeReturnType | transpose () const |
|
const AdjointReturnType | adjoint () const |
|
constexpr EIGEN_DEVICE_FUNC Derived & | derived () |
|
constexpr EIGEN_DEVICE_FUNC const Derived & | derived () const |
|
constexpr EIGEN_DEVICE_FUNC Derived & | derived () |
|
constexpr EIGEN_DEVICE_FUNC const Derived & | derived () const |
|
EIGEN_DEVICE_FUNC Derived & | const_cast_derived () const |
|
EIGEN_DEVICE_FUNC const Derived & | const_derived () const |
|
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index | rows () const EIGEN_NOEXCEPT |
|
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index | cols () const EIGEN_NOEXCEPT |
|
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index | size () const EIGEN_NOEXCEPT |
|
template<typename Dest > |
EIGEN_DEVICE_FUNC void | evalTo (Dest &dst) const |
|
template<typename Dest > |
EIGEN_DEVICE_FUNC void | addTo (Dest &dst) const |
|
template<typename Dest > |
EIGEN_DEVICE_FUNC void | subTo (Dest &dst) const |
|
template<typename Dest > |
EIGEN_DEVICE_FUNC void | applyThisOnTheRight (Dest &dst) const |
|
template<typename Dest > |
EIGEN_DEVICE_FUNC void | applyThisOnTheLeft (Dest &dst) const |
|
template<typename Device > |
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE DeviceWrapper< Derived, Device > | device (Device &device) |
|
template<typename Device > |
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE DeviceWrapper< const Derived, Device > | device (Device &device) const |
|
template<typename MatrixType_, int Options_>
class Eigen::JacobiSVD< MatrixType_, Options_ >
Two-sided Jacobi SVD decomposition of a rectangular matrix.
- Template Parameters
-
MatrixType_ | the type of the matrix of which we are computing the SVD decomposition |
Options | this optional parameter allows one to specify the type of QR decomposition that will be used internally for the R-SVD step for non-square matrices. Additionally, it allows one to specify whether to compute thin or full unitaries U and V. See discussion of possible values below. |
SVD decomposition consists in decomposing any n-by-p matrix A as a product
\[ A = U S V^* \]
where U is a n-by-n unitary, V is a p-by-p unitary, and S is a n-by-p real positive matrix which is zero outside of its main diagonal; the diagonal entries of S are known as the singular values of A and the columns of U and V are known as the left and right singular vectors of A respectively.
Singular values are always sorted in decreasing order.
This JacobiSVD decomposition computes only the singular values by default. If you want U or V, you need to ask for them explicitly.
You can ask for only thin U or V to be computed, meaning the following. In case of a rectangular n-by-p matrix, letting m be the smaller value among n and p, there are only m singular vectors; the remaining columns of U and V do not correspond to actual singular vectors. Asking for thin U or V means asking for only their m first columns to be formed. So U is then a n-by-m matrix, and V is then a p-by-m matrix. Notice that thin U and V are all you need for (least squares) solving.
Here's an example demonstrating basic usage:
MatrixXf
m = MatrixXf::Random(3, 2);
cout <<
"Here is the matrix m:" << endl <<
m << endl;
JacobiSVD<MatrixXf, ComputeThinU | ComputeThinV>
svd(
m);
cout <<
"Its singular values are:" << endl <<
svd.singularValues() << endl;
cout <<
"Its left singular vectors are the columns of the thin U matrix:" << endl <<
svd.matrixU() << endl;
cout <<
"Its right singular vectors are the columns of the thin V matrix:" << endl <<
svd.matrixV() << endl;
Vector3f rhs(1, 0, 0);
cout << "Now consider this rhs vector:" << endl << rhs << endl;
cout <<
"A least-squares solution of m*x = rhs is:" << endl <<
svd.solve(rhs) << endl;
cout<< "Here is the matrix m:"<< endl<< m<< endl;JacobiSVD< MatrixXf, ComputeThinU|ComputeThinV > svd(m)
int * m
Definition: level2_cplx_impl.h:294
Output:
This JacobiSVD class is a two-sided Jacobi R-SVD decomposition, ensuring optimal reliability and accuracy. The downside is that it's slower than bidiagonalizing SVD algorithms for large square matrices; however its complexity is still \( O(n^2p) \) where n is the smaller dimension and p is the greater dimension, meaning that it is still of the same order of complexity as the faster bidiagonalizing R-SVD algorithms. In particular, like any R-SVD, it takes advantage of non-squareness in that its complexity is only linear in the greater dimension.
If the input matrix has inf or nan coefficients, the result of the computation is undefined, but the computation is guaranteed to terminate in finite (and reasonable) time.
The possible QR preconditioners that can be set with Options template parameter are:
- ColPivHouseholderQRPreconditioner is the default. In practice it's very safe. It uses column-pivoting QR.
- FullPivHouseholderQRPreconditioner, is the safest and slowest. It uses full-pivoting QR. Contrary to other QRs, it doesn't allow computing thin unitaries.
- HouseholderQRPreconditioner is the fastest, and less safe and accurate than the pivoting variants. It uses non-pivoting QR. This is very similar in safety and accuracy to the bidiagonalization process used by bidiagonalizing SVD algorithms (since bidiagonalization is inherently non-pivoting). However the resulting SVD is still more reliable than bidiagonalizing SVDs because the Jacobi-based iterarive process is more reliable than the optimized bidiagonal SVD iterations.
- NoQRPreconditioner allows not to use a QR preconditioner at all. This is useful if you know that you will only be computing JacobiSVD decompositions of square matrices. Non-square matrices require a QR preconditioner. Using this option will result in faster compilation and smaller executable code. It won't significantly speed up computation, since JacobiSVD is always checking if QR preconditioning is needed before applying it anyway.
One may also use the Options template parameter to specify how the unitaries should be computed. The options are ComputeThinU, ComputeThinV, ComputeFullU, ComputeFullV. It is not possible to request both the thin and full versions of a unitary. By default, unitaries will not be computed.
You can set the QRPreconditioner and unitary options together: JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner | ComputeThinU | ComputeFullV>
- See also
- MatrixBase::jacobiSvd()
template<typename MatrixType , int Options>
695 matrix.template topLeftCorner<DiagSizeAtCompileTime, DiagSizeAtCompileTime>(
diagSize(),
diagSize()) / scale;
723 JacobiRotation<RealScalar> j_left, j_right;
734 maxDiagEntry = numext::maxi<RealScalar>(
AnnoyingScalar abs(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:135
AnnoyingScalar imag(const AnnoyingScalar &)
Definition: AnnoyingScalar.h:132
int i
Definition: BiCGSTAB_step_by_step.cpp:9
float * p
Definition: Tutorial_Map_using.cpp:9
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:46
Base::RealScalar RealScalar
Definition: JacobiSVD.h:506
EIGEN_DEVICE_FUNC constexpr EIGEN_STRONG_INLINE const Scalar & coeff(Index rowId, Index colId) const
Definition: PlainObjectBase.h:198
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void swap(DenseBase< OtherDerived > &other)
Override DenseBase::swap() since for dynamic-sized matrices of same type it is enough to swap the dat...
Definition: PlainObjectBase.h:898
bool m_computeFullV
Definition: SVDBase.h:337
ComputationInfo m_info
Definition: SVDBase.h:334
bool m_computeThinU
Definition: SVDBase.h:336
bool computeV() const
Definition: SVDBase.h:275
bool m_isInitialized
Definition: SVDBase.h:335
MatrixVType m_matrixV
Definition: SVDBase.h:332
bool computeU() const
Definition: SVDBase.h:273
RealScalar threshold() const
Definition: SVDBase.h:265
Index m_nonzeroSingularValues
Definition: SVDBase.h:339
SingularValuesType m_singularValues
Definition: SVDBase.h:333
bool m_computeThinV
Definition: SVDBase.h:337
bool m_computeFullU
Definition: SVDBase.h:336
MatrixUType m_matrixU
Definition: SVDBase.h:331
float real
Definition: datatypes.h:10
#define min(a, b)
Definition: datatypes.h:22
@ InvalidInput
Definition: Constants.h:447
EIGEN_BLAS_FUNC() swap(int *n, RealScalar *px, int *incx, RealScalar *py, int *incy)
Definition: level1_impl.h:117
const Scalar * a
Definition: level2_cplx_impl.h:32
void real_2x2_jacobi_svd(const MatrixType &matrix, Index p, Index q, JacobiRotation< RealScalar > *j_left, JacobiRotation< RealScalar > *j_right)
Definition: RealSvd2x2.h:22
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE bool() isfinite(const Eigen::bfloat16 &h)
Definition: BFloat16.h:752
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bool is_exactly_zero(const X &x)
Definition: Meta.h:592
EIGEN_DEVICE_FUNC const Scalar & q
Definition: SpecialFunctionsImpl.h:2019
bool finished
Definition: MergeRestartFiles.py:79
double epsilon
Definition: osc_ring_sarah_asymptotics.h:43
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:43
@ IsComplex
Definition: NumTraits.h:176
void run(const string &dir_name, LinearSolver *linear_solver_pt, const unsigned nel_1d, bool mess_up_order)
Definition: two_d_poisson_compare_solvers.cc:317
References a, abs(), Eigen::JacobiSVD< MatrixType_, Options_ >::allocate(), Eigen::PlainObjectBase< Derived >::coeff(), Eigen::JacobiSVD< MatrixType_, Options_ >::cols(), Eigen::SVDBase< JacobiSVD< MatrixType_, Options_ > >::computeU(), Eigen::SVDBase< JacobiSVD< MatrixType_, Options_ > >::computeV(), Eigen::SVDBase< JacobiSVD< MatrixType_, Options_ > >::diagSize(), oomph::SarahBL::epsilon, MergeRestartFiles::finished, i, imag(), Eigen::InvalidInput, Eigen::numext::is_exactly_zero(), Eigen::numext::isfinite(), Eigen::SVDBase< JacobiSVD< MatrixType_, Options_ > >::m_computeFullU, Eigen::SVDBase< JacobiSVD< MatrixType_, Options_ > >::m_computeFullV, Eigen::SVDBase< JacobiSVD< MatrixType_, Options_ > >::m_computeThinU, Eigen::SVDBase< JacobiSVD< MatrixType_, Options_ > >::m_computeThinV, Eigen::SVDBase< JacobiSVD< MatrixType_, Options_ > >::m_info, Eigen::SVDBase< JacobiSVD< MatrixType_, Options_ > >::m_isInitialized, Eigen::SVDBase< JacobiSVD< MatrixType_, Options_ > >::m_matrixU, Eigen::SVDBase< JacobiSVD< MatrixType_, Options_ > >::m_matrixV, Eigen::SVDBase< JacobiSVD< MatrixType_, Options_ > >::m_nonzeroSingularValues, Eigen::JacobiSVD< MatrixType_, Options_ >::m_qr_precond_morecols, Eigen::JacobiSVD< MatrixType_, Options_ >::m_qr_precond_morerows, Eigen::SVDBase< JacobiSVD< MatrixType_, Options_ > >::m_singularValues, Eigen::JacobiSVD< MatrixType_, Options_ >::m_workMatrix, matrix(), min, p, Eigen::numext::q, Eigen::internal::real_2x2_jacobi_svd(), Eigen::JacobiSVD< MatrixType_, Options_ >::rows(), Eigen::PlainObjectBase< Derived >::swap(), swap(), Eigen::SVDBase< JacobiSVD< MatrixType_, Options_ > >::threshold(), and Eigen::JacobiRotation< Scalar >::transpose().
Referenced by Eigen::JacobiSVD< MatrixType_, Options_ >::compute(), and Eigen::JacobiSVD< MatrixType_, Options_ >::JacobiSVD().