MOAB
4.9.3pre
|
Two-sided Jacobi SVD decomposition of a rectangular matrix. More...
#include <JacobiSVD.h>
Public Types | |
enum | { RowsAtCompileTime = MatrixType::RowsAtCompileTime, ColsAtCompileTime = MatrixType::ColsAtCompileTime, DiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_DYNAMIC(RowsAtCompileTime,ColsAtCompileTime), MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime, MaxDiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(MaxRowsAtCompileTime,MaxColsAtCompileTime), MatrixOptions = MatrixType::Options } |
typedef _MatrixType | MatrixType |
typedef MatrixType::Scalar | Scalar |
typedef NumTraits< typename MatrixType::Scalar >::Real | RealScalar |
typedef Base::MatrixUType | MatrixUType |
typedef Base::MatrixVType | MatrixVType |
typedef Base::SingularValuesType | SingularValuesType |
typedef internal::plain_row_type < MatrixType >::type | RowType |
typedef internal::plain_col_type < MatrixType >::type | ColType |
typedef Matrix< Scalar, DiagSizeAtCompileTime, DiagSizeAtCompileTime, MatrixOptions, MaxDiagSizeAtCompileTime, MaxDiagSizeAtCompileTime > | WorkMatrixType |
Public Member Functions | |
JacobiSVD () | |
Default Constructor. | |
JacobiSVD (Index p_rows, Index p_cols, unsigned int computationOptions=0) | |
Default Constructor with memory preallocation. | |
JacobiSVD (const MatrixType &matrix, unsigned int computationOptions=0) | |
Constructor performing the decomposition of given matrix. | |
JacobiSVD & | compute (const MatrixType &matrix, unsigned int computationOptions) |
Method performing the decomposition of given matrix using custom options. | |
JacobiSVD & | compute (const MatrixType &matrix) |
Method performing the decomposition of given matrix using current options. | |
Protected Attributes | |
WorkMatrixType | m_workMatrix |
internal::qr_preconditioner_impl < MatrixType, QRPreconditioner, internal::PreconditionIfMoreColsThanRows > | m_qr_precond_morecols |
internal::qr_preconditioner_impl < MatrixType, QRPreconditioner, internal::PreconditionIfMoreRowsThanCols > | m_qr_precond_morerows |
MatrixType | m_scaledMatrix |
Private Types | |
typedef SVDBase< JacobiSVD > | Base |
Private Member Functions | |
void | allocate (Index rows, Index cols, unsigned int computationOptions) |
Friends | |
struct | internal::svd_precondition_2x2_block_to_be_real |
struct | internal::qr_preconditioner_impl |
Two-sided Jacobi SVD decomposition of a rectangular matrix.
_MatrixType | the type of the matrix of which we are computing the SVD decomposition |
QRPreconditioner | this optional parameter allows to specify the type of QR decomposition that will be used internally for the R-SVD step for non-square matrices. See discussion of possible values below. |
SVD decomposition consists in decomposing any n-by-p matrix A as a product
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:
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 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 values for QRPreconditioner are:
Definition at line 498 of file JacobiSVD.h.
typedef SVDBase<JacobiSVD> Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::Base [private] |
Definition at line 501 of file JacobiSVD.h.
typedef internal::plain_col_type<MatrixType>::type Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::ColType |
Definition at line 522 of file JacobiSVD.h.
typedef _MatrixType Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::MatrixType |
Reimplemented from Eigen::SVDBase< JacobiSVD< _MatrixType, QRPreconditioner > >.
Definition at line 504 of file JacobiSVD.h.
typedef Base::MatrixUType Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::MatrixUType |
Reimplemented from Eigen::SVDBase< JacobiSVD< _MatrixType, QRPreconditioner > >.
Definition at line 517 of file JacobiSVD.h.
typedef Base::MatrixVType Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::MatrixVType |
Reimplemented from Eigen::SVDBase< JacobiSVD< _MatrixType, QRPreconditioner > >.
Definition at line 518 of file JacobiSVD.h.
typedef NumTraits<typename MatrixType::Scalar>::Real Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::RealScalar |
Reimplemented from Eigen::SVDBase< JacobiSVD< _MatrixType, QRPreconditioner > >.
Definition at line 506 of file JacobiSVD.h.
typedef internal::plain_row_type<MatrixType>::type Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::RowType |
Definition at line 521 of file JacobiSVD.h.
typedef MatrixType::Scalar Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::Scalar |
Reimplemented from Eigen::SVDBase< JacobiSVD< _MatrixType, QRPreconditioner > >.
Definition at line 505 of file JacobiSVD.h.
typedef Base::SingularValuesType Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::SingularValuesType |
Reimplemented from Eigen::SVDBase< JacobiSVD< _MatrixType, QRPreconditioner > >.
Definition at line 519 of file JacobiSVD.h.
typedef Matrix<Scalar, DiagSizeAtCompileTime, DiagSizeAtCompileTime, MatrixOptions, MaxDiagSizeAtCompileTime, MaxDiagSizeAtCompileTime> Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::WorkMatrixType |
Definition at line 525 of file JacobiSVD.h.
anonymous enum |
RowsAtCompileTime | |
ColsAtCompileTime | |
DiagSizeAtCompileTime | |
MaxRowsAtCompileTime | |
MaxColsAtCompileTime | |
MaxDiagSizeAtCompileTime | |
MatrixOptions |
Definition at line 507 of file JacobiSVD.h.
{ RowsAtCompileTime = MatrixType::RowsAtCompileTime, ColsAtCompileTime = MatrixType::ColsAtCompileTime, DiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_DYNAMIC(RowsAtCompileTime,ColsAtCompileTime), MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime, MaxDiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(MaxRowsAtCompileTime,MaxColsAtCompileTime), MatrixOptions = MatrixType::Options };
Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::JacobiSVD | ( | ) | [inline] |
Default Constructor.
The default constructor is useful in cases in which the user intends to perform decompositions via JacobiSVD::compute(const MatrixType&).
Definition at line 532 of file JacobiSVD.h.
{}
Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::JacobiSVD | ( | Index | p_rows, |
Index | p_cols, | ||
unsigned int | computationOptions = 0 |
||
) | [inline] |
Default Constructor with memory preallocation.
Like the default constructor but with preallocation of the internal data according to the specified problem size.
Definition at line 542 of file JacobiSVD.h.
{ allocate(p_rows, p_cols, computationOptions); }
Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::JacobiSVD | ( | const MatrixType & | matrix, |
unsigned int | computationOptions = 0 |
||
) | [inline, explicit] |
Constructor performing the decomposition of given matrix.
matrix | the matrix to decompose |
computationOptions | optional parameter allowing to specify if you want full or thin U or V unitaries to be computed. By default, none is computed. This is a bit-field, the possible bits are ComputeFullU, ComputeThinU, ComputeFullV, ComputeThinV. |
Thin unitaries are only available if your matrix type has a Dynamic number of columns (for example MatrixXf). They also are not available with the (non-default) FullPivHouseholderQR preconditioner.
Definition at line 557 of file JacobiSVD.h.
{ compute(matrix, computationOptions); }
void Eigen::JacobiSVD< MatrixType, QRPreconditioner >::allocate | ( | Index | rows, |
Index | cols, | ||
unsigned int | computationOptions | ||
) | [private] |
Reimplemented from Eigen::SVDBase< JacobiSVD< _MatrixType, QRPreconditioner > >.
Definition at line 624 of file JacobiSVD.h.
{ eigen_assert(p_rows >= 0 && p_cols >= 0); if (m_isAllocated && p_rows == m_rows && p_cols == m_cols && computationOptions == m_computationOptions) { return; } m_rows = p_rows; m_cols = p_cols; m_isInitialized = false; m_isAllocated = true; m_computationOptions = computationOptions; m_computeFullU = (computationOptions & ComputeFullU) != 0; m_computeThinU = (computationOptions & ComputeThinU) != 0; m_computeFullV = (computationOptions & ComputeFullV) != 0; m_computeThinV = (computationOptions & ComputeThinV) != 0; eigen_assert(!(m_computeFullU && m_computeThinU) && "JacobiSVD: you can't ask for both full and thin U"); eigen_assert(!(m_computeFullV && m_computeThinV) && "JacobiSVD: you can't ask for both full and thin V"); eigen_assert(EIGEN_IMPLIES(m_computeThinU || m_computeThinV, MatrixType::ColsAtCompileTime==Dynamic) && "JacobiSVD: thin U and V are only available when your matrix has a dynamic number of columns."); if (QRPreconditioner == FullPivHouseholderQRPreconditioner) { eigen_assert(!(m_computeThinU || m_computeThinV) && "JacobiSVD: can't compute thin U or thin V with the FullPivHouseholderQR preconditioner. " "Use the ColPivHouseholderQR preconditioner instead."); } m_diagSize = (std::min)(m_rows, m_cols); m_singularValues.resize(m_diagSize); if(RowsAtCompileTime==Dynamic) m_matrixU.resize(m_rows, m_computeFullU ? m_rows : m_computeThinU ? m_diagSize : 0); if(ColsAtCompileTime==Dynamic) m_matrixV.resize(m_cols, m_computeFullV ? m_cols : m_computeThinV ? m_diagSize : 0); m_workMatrix.resize(m_diagSize, m_diagSize); if(m_cols>m_rows) m_qr_precond_morecols.allocate(*this); if(m_rows>m_cols) m_qr_precond_morerows.allocate(*this); if(m_rows!=m_cols) m_scaledMatrix.resize(p_rows,p_cols); }
JacobiSVD< MatrixType, QRPreconditioner > & Eigen::JacobiSVD< MatrixType, QRPreconditioner >::compute | ( | const MatrixType & | matrix, |
unsigned int | computationOptions | ||
) |
Method performing the decomposition of given matrix using custom options.
matrix | the matrix to decompose |
computationOptions | optional parameter allowing to specify if you want full or thin U or V unitaries to be computed. By default, none is computed. This is a bit-field, the possible bits are ComputeFullU, ComputeThinU, ComputeFullV, ComputeThinV. |
Thin unitaries are only available if your matrix type has a Dynamic number of columns (for example MatrixXf). They also are not available with the (non-default) FullPivHouseholderQR preconditioner.
Definition at line 674 of file JacobiSVD.h.
{ using std::abs; allocate(matrix.rows(), matrix.cols(), computationOptions); // currently we stop when we reach precision 2*epsilon as the last bit of precision can require an unreasonable number of iterations, // only worsening the precision of U and V as we accumulate more rotations const RealScalar precision = RealScalar(2) * NumTraits<Scalar>::epsilon(); // limit for very small denormal numbers to be considered zero in order to avoid infinite loops (see bug 286) // FIXME What about considerering any denormal numbers as zero, using: // const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)(); const RealScalar considerAsZero = RealScalar(2) * std::numeric_limits<RealScalar>::denorm_min(); // Scaling factor to reduce over/under-flows RealScalar scale = matrix.cwiseAbs().maxCoeff(); if(scale==RealScalar(0)) scale = RealScalar(1); /*** step 1. The R-SVD step: we use a QR decomposition to reduce to the case of a square matrix */ if(m_rows!=m_cols) { m_scaledMatrix = matrix / scale; m_qr_precond_morecols.run(*this, m_scaledMatrix); m_qr_precond_morerows.run(*this, m_scaledMatrix); } else { m_workMatrix = matrix.block(0,0,m_diagSize,m_diagSize) / scale; if(m_computeFullU) m_matrixU.setIdentity(m_rows,m_rows); if(m_computeThinU) m_matrixU.setIdentity(m_rows,m_diagSize); if(m_computeFullV) m_matrixV.setIdentity(m_cols,m_cols); if(m_computeThinV) m_matrixV.setIdentity(m_cols, m_diagSize); } /*** step 2. The main Jacobi SVD iteration. ***/ bool finished = false; while(!finished) { finished = true; // do a sweep: for all index pairs (p,q), perform SVD of the corresponding 2x2 sub-matrix for(Index p = 1; p < m_diagSize; ++p) { for(Index q = 0; q < p; ++q) { // if this 2x2 sub-matrix is not diagonal already... // notice that this comparison will evaluate to false if any NaN is involved, ensuring that NaN's don't // keep us iterating forever. Similarly, small denormal numbers are considered zero. RealScalar threshold = numext::maxi<RealScalar>(considerAsZero, precision * numext::maxi<RealScalar>(abs(m_workMatrix.coeff(p,p)), abs(m_workMatrix.coeff(q,q)))); // We compare both values to threshold instead of calling max to be robust to NaN (See bug 791) if(abs(m_workMatrix.coeff(p,q))>threshold || abs(m_workMatrix.coeff(q,p)) > threshold) { finished = false; // perform SVD decomposition of 2x2 sub-matrix corresponding to indices p,q to make it diagonal internal::svd_precondition_2x2_block_to_be_real<MatrixType, QRPreconditioner>::run(m_workMatrix, *this, p, q); JacobiRotation<RealScalar> j_left, j_right; internal::real_2x2_jacobi_svd(m_workMatrix, p, q, &j_left, &j_right); // accumulate resulting Jacobi rotations m_workMatrix.applyOnTheLeft(p,q,j_left); if(computeU()) m_matrixU.applyOnTheRight(p,q,j_left.transpose()); m_workMatrix.applyOnTheRight(p,q,j_right); if(computeV()) m_matrixV.applyOnTheRight(p,q,j_right); } } } } /*** step 3. The work matrix is now diagonal, so ensure it's positive so its diagonal entries are the singular values ***/ for(Index i = 0; i < m_diagSize; ++i) { RealScalar a = abs(m_workMatrix.coeff(i,i)); m_singularValues.coeffRef(i) = a; if(computeU() && (a!=RealScalar(0))) m_matrixU.col(i) *= m_workMatrix.coeff(i,i)/a; } m_singularValues *= scale; /*** step 4. Sort singular values in descending order and compute the number of nonzero singular values ***/ m_nonzeroSingularValues = m_diagSize; for(Index i = 0; i < m_diagSize; i++) { Index pos; RealScalar maxRemainingSingularValue = m_singularValues.tail(m_diagSize-i).maxCoeff(&pos); if(maxRemainingSingularValue == RealScalar(0)) { m_nonzeroSingularValues = i; break; } if(pos) { pos += i; std::swap(m_singularValues.coeffRef(i), m_singularValues.coeffRef(pos)); if(computeU()) m_matrixU.col(pos).swap(m_matrixU.col(i)); if(computeV()) m_matrixV.col(pos).swap(m_matrixV.col(i)); } } m_isInitialized = true; return *this; }
JacobiSVD& Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::compute | ( | const MatrixType & | matrix | ) | [inline] |
Method performing the decomposition of given matrix using current options.
matrix | the matrix to decompose |
This method uses the current computationOptions, as already passed to the constructor or to compute(const MatrixType&, unsigned int).
Definition at line 580 of file JacobiSVD.h.
{ return compute(matrix, m_computationOptions); }
friend struct internal::qr_preconditioner_impl [friend] |
Definition at line 616 of file JacobiSVD.h.
friend struct internal::svd_precondition_2x2_block_to_be_real [friend] |
Definition at line 614 of file JacobiSVD.h.
internal::qr_preconditioner_impl<MatrixType, QRPreconditioner, internal::PreconditionIfMoreColsThanRows> Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::m_qr_precond_morecols [protected] |
Definition at line 618 of file JacobiSVD.h.
internal::qr_preconditioner_impl<MatrixType, QRPreconditioner, internal::PreconditionIfMoreRowsThanCols> Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::m_qr_precond_morerows [protected] |
Definition at line 619 of file JacobiSVD.h.
MatrixType Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::m_scaledMatrix [protected] |
Definition at line 620 of file JacobiSVD.h.
WorkMatrixType Eigen::JacobiSVD< _MatrixType, QRPreconditioner >::m_workMatrix [protected] |
Definition at line 611 of file JacobiSVD.h.