10 #ifndef EIGEN_JACOBISVD_H 11 #define EIGEN_JACOBISVD_H 18 template<
typename MatrixType,
int QRPreconditioner,
31 template<
typename MatrixType,
int QRPreconditioner,
int Case>
34 enum { a = MatrixType::RowsAtCompileTime !=
Dynamic &&
35 MatrixType::ColsAtCompileTime !=
Dynamic &&
36 MatrixType::ColsAtCompileTime <= MatrixType::RowsAtCompileTime,
37 b = MatrixType::RowsAtCompileTime !=
Dynamic &&
38 MatrixType::ColsAtCompileTime !=
Dynamic &&
39 MatrixType::RowsAtCompileTime <= MatrixType::ColsAtCompileTime,
42 (Case == PreconditionIfMoreRowsThanCols &&
bool(b)) )
46 template<
typename MatrixType,
int QRPreconditioner,
int Case,
50 template<
typename MatrixType,
int QRPreconditioner,
int Case>
54 typedef typename MatrixType::Index
Index;
64 template<
typename MatrixType>
68 typedef typename MatrixType::Index
Index;
69 typedef typename MatrixType::Scalar
Scalar;
72 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
73 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime
79 if (svd.
rows() != m_qr.rows() || svd.
cols() != m_qr.cols())
89 if(matrix.rows() > matrix.cols())
92 svd.
m_workMatrix = m_qr.matrixQR().block(0,0,matrix.cols(),matrix.cols()).
template triangularView<Upper>();
105 template<
typename MatrixType>
109 typedef typename MatrixType::Index
Index;
110 typedef typename MatrixType::Scalar
Scalar;
113 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
114 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
115 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
116 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
117 Options = MatrixType::Options
124 if (svd.
cols() != m_qr.rows() || svd.
rows() != m_qr.cols())
129 m_adjoint.resize(svd.
cols(), svd.
rows());
135 if(matrix.cols() > matrix.rows())
137 m_adjoint = matrix.adjoint();
138 m_qr.compute(m_adjoint);
139 svd.
m_workMatrix = m_qr.matrixQR().block(0,0,matrix.rows(),matrix.rows()).
template triangularView<Upper>().adjoint();
155 template<
typename MatrixType>
159 typedef typename MatrixType::Index
Index;
163 if (svd.
rows() != m_qr.rows() || svd.
cols() != m_qr.cols())
174 if(matrix.rows() > matrix.cols())
176 m_qr.compute(matrix);
177 svd.
m_workMatrix = m_qr.matrixQR().block(0,0,matrix.cols(),matrix.cols()).
template triangularView<Upper>();
181 svd.
m_matrixU.setIdentity(matrix.rows(), matrix.cols());
182 m_qr.householderQ().applyThisOnTheLeft(svd.
m_matrixU, m_workspace);
196 template<
typename MatrixType>
200 typedef typename MatrixType::Index
Index;
201 typedef typename MatrixType::Scalar
Scalar;
204 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
205 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
206 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
207 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
208 Options = MatrixType::Options
216 if (svd.
cols() != m_qr.rows() || svd.
rows() != m_qr.cols())
223 m_adjoint.resize(svd.
cols(), svd.
rows());
228 if(matrix.cols() > matrix.rows())
230 m_adjoint = matrix.adjoint();
231 m_qr.compute(m_adjoint);
233 svd.
m_workMatrix = m_qr.matrixQR().block(0,0,matrix.rows(),matrix.rows()).
template triangularView<Upper>().adjoint();
237 svd.
m_matrixV.setIdentity(matrix.cols(), matrix.rows());
238 m_qr.householderQ().applyThisOnTheLeft(svd.
m_matrixV, m_workspace);
255 template<
typename MatrixType>
259 typedef typename MatrixType::Index
Index;
263 if (svd.
rows() != m_qr.rows() || svd.
cols() != m_qr.cols())
274 if(matrix.rows() > matrix.cols())
276 m_qr.compute(matrix);
277 svd.
m_workMatrix = m_qr.matrixQR().block(0,0,matrix.cols(),matrix.cols()).
template triangularView<Upper>();
281 svd.
m_matrixU.setIdentity(matrix.rows(), matrix.cols());
282 m_qr.householderQ().applyThisOnTheLeft(svd.
m_matrixU, m_workspace);
295 template<
typename MatrixType>
299 typedef typename MatrixType::Index
Index;
300 typedef typename MatrixType::Scalar
Scalar;
303 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
304 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
305 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
306 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
307 Options = MatrixType::Options
315 if (svd.
cols() != m_qr.rows() || svd.
rows() != m_qr.cols())
322 m_adjoint.resize(svd.
cols(), svd.
rows());
327 if(matrix.cols() > matrix.rows())
329 m_adjoint = matrix.adjoint();
330 m_qr.compute(m_adjoint);
332 svd.
m_workMatrix = m_qr.matrixQR().block(0,0,matrix.rows(),matrix.rows()).
template triangularView<Upper>().adjoint();
336 svd.
m_matrixV.setIdentity(matrix.cols(), matrix.rows());
337 m_qr.householderQ().applyThisOnTheLeft(svd.
m_matrixV, m_workspace);
357 template<
typename MatrixType,
int QRPreconditioner>
365 template<
typename MatrixType,
int QRPreconditioner>
369 typedef typename MatrixType::Scalar
Scalar;
381 work_matrix.row(p) *= z;
384 work_matrix.row(q) *= z;
389 rot.
c() = conj(work_matrix.
coeff(p,p)) / n;
390 rot.
s() = work_matrix.
coeff(q,p) / n;
391 work_matrix.applyOnTheLeft(p,q,rot);
393 if(work_matrix.
coeff(p,q) != Scalar(0))
395 Scalar z =
abs(work_matrix.
coeff(p,q)) / work_matrix.
coeff(p,q);
396 work_matrix.col(q) *= z;
399 if(work_matrix.
coeff(q,q) != Scalar(0))
402 work_matrix.row(q) *= z;
409 template<
typename MatrixType,
typename RealScalar,
typename Index>
421 if(t == RealScalar(0))
423 rot1.c() = RealScalar(0);
424 rot1.s() = d > RealScalar(0) ? RealScalar(1) : RealScalar(-1);
428 RealScalar u = d / t;
430 rot1.s() = rot1.c() * u;
432 m.applyOnTheLeft(0,1,rot1);
492 template<
typename _MatrixType,
int QRPreconditioner>
class JacobiSVD 497 typedef typename MatrixType::Scalar
Scalar;
499 typedef typename MatrixType::Index
Index;
501 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
502 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
504 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
505 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
507 MatrixOptions = MatrixType::Options
510 typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime,
511 MatrixOptions, MaxRowsAtCompileTime, MaxRowsAtCompileTime>
513 typedef Matrix<Scalar, ColsAtCompileTime, ColsAtCompileTime,
514 MatrixOptions, MaxColsAtCompileTime, MaxColsAtCompileTime>
519 typedef Matrix<Scalar, DiagSizeAtCompileTime, DiagSizeAtCompileTime,
520 MatrixOptions, MaxDiagSizeAtCompileTime, MaxDiagSizeAtCompileTime>
529 : m_isInitialized(false),
530 m_isAllocated(false),
531 m_computationOptions(0),
532 m_rows(-1), m_cols(-1)
542 JacobiSVD(Index rows, Index cols,
unsigned int computationOptions = 0)
543 : m_isInitialized(false),
544 m_isAllocated(false),
545 m_computationOptions(0),
546 m_rows(-1), m_cols(-1)
548 allocate(rows, cols, computationOptions);
561 JacobiSVD(
const MatrixType& matrix,
unsigned int computationOptions = 0)
562 : m_isInitialized(false),
563 m_isAllocated(false),
564 m_computationOptions(0),
565 m_rows(-1), m_cols(-1)
567 compute(matrix, computationOptions);
580 JacobiSVD& compute(
const MatrixType& matrix,
unsigned int computationOptions);
590 return compute(matrix, m_computationOptions);
604 eigen_assert(m_isInitialized &&
"JacobiSVD is not initialized.");
605 eigen_assert(computeU() &&
"This JacobiSVD decomposition didn't compute U. Did you ask for it?");
620 eigen_assert(m_isInitialized &&
"JacobiSVD is not initialized.");
621 eigen_assert(computeV() &&
"This JacobiSVD decomposition didn't compute V. Did you ask for it?");
632 eigen_assert(m_isInitialized &&
"JacobiSVD is not initialized.");
633 return m_singularValues;
637 inline bool computeU()
const {
return m_computeFullU || m_computeThinU; }
639 inline bool computeV()
const {
return m_computeFullV || m_computeThinV; }
650 template<
typename Rhs>
654 eigen_assert(m_isInitialized &&
"JacobiSVD is not initialized.");
655 eigen_assert(computeU() && computeV() &&
"JacobiSVD::solve() requires both unitaries U and V to be computed (thin unitaries suffice).");
662 eigen_assert(m_isInitialized &&
"JacobiSVD is not initialized.");
663 return m_nonzeroSingularValues;
666 inline Index
rows()
const {
return m_rows; }
667 inline Index
cols()
const {
return m_cols; }
670 void allocate(Index rows, Index cols,
unsigned int computationOptions);
681 Index m_nonzeroSingularValues,
m_rows, m_cols, m_diagSize;
683 template<
typename __MatrixType,
int _QRPreconditioner,
bool _IsComplex>
685 template<
typename __MatrixType,
int _QRPreconditioner,
int _Case,
bool _DoAnything>
692 template<
typename MatrixType,
int QRPreconditioner>
700 computationOptions == m_computationOptions)
707 m_isInitialized =
false;
708 m_isAllocated =
true;
709 m_computationOptions = computationOptions;
710 m_computeFullU = (computationOptions &
ComputeFullU) != 0;
711 m_computeThinU = (computationOptions &
ComputeThinU) != 0;
712 m_computeFullV = (computationOptions &
ComputeFullV) != 0;
713 m_computeThinV = (computationOptions &
ComputeThinV) != 0;
714 eigen_assert(!(m_computeFullU && m_computeThinU) &&
"JacobiSVD: you can't ask for both full and thin U");
715 eigen_assert(!(m_computeFullV && m_computeThinV) &&
"JacobiSVD: you can't ask for both full and thin V");
717 "JacobiSVD: thin U and V are only available when your matrix has a dynamic number of columns.");
721 "JacobiSVD: can't compute thin U or thin V with the FullPivHouseholderQR preconditioner. " 722 "Use the ColPivHouseholderQR preconditioner instead.");
724 m_diagSize = (std::min)(m_rows, m_cols);
725 m_singularValues.resize(m_diagSize);
727 m_matrixU.resize(m_rows, m_computeFullU ? m_rows
728 : m_computeThinU ? m_diagSize
731 m_matrixV.resize(m_cols, m_computeFullV ? m_cols
732 : m_computeThinV ? m_diagSize
734 m_workMatrix.resize(m_diagSize, m_diagSize);
736 if(m_cols>m_rows) m_qr_precond_morecols.allocate(*
this);
737 if(m_rows>m_cols) m_qr_precond_morerows.allocate(*
this);
740 template<
typename MatrixType,
int QRPreconditioner>
745 allocate(matrix.rows(), matrix.cols(), computationOptions);
752 const RealScalar considerAsZero =
RealScalar(2) * std::numeric_limits<RealScalar>::denorm_min();
756 if(!m_qr_precond_morecols.run(*
this, matrix) && !m_qr_precond_morerows.run(*
this, matrix))
758 m_workMatrix = matrix.block(0,0,m_diagSize,m_diagSize);
759 if(m_computeFullU) m_matrixU.setIdentity(m_rows,m_rows);
760 if(m_computeThinU) m_matrixU.setIdentity(m_rows,m_diagSize);
761 if(m_computeFullV) m_matrixV.setIdentity(m_cols,m_cols);
762 if(m_computeThinV) m_matrixV.setIdentity(m_cols, m_diagSize);
767 bool finished =
false;
774 for(
Index p = 1; p < m_diagSize; ++p)
776 for(
Index q = 0; q < p; ++q)
782 RealScalar threshold = (max)(considerAsZero, precision * (max)(
abs(m_workMatrix.coeff(p,p)),
783 abs(m_workMatrix.coeff(q,q))));
784 if((max)(
abs(m_workMatrix.coeff(p,q)),
abs(m_workMatrix.coeff(q,p))) > threshold)
794 m_workMatrix.applyOnTheLeft(p,q,j_left);
795 if(computeU()) m_matrixU.applyOnTheRight(p,q,j_left.
transpose());
797 m_workMatrix.applyOnTheRight(p,q,j_right);
798 if(computeV()) m_matrixV.applyOnTheRight(p,q,j_right);
806 for(
Index i = 0; i < m_diagSize; ++i)
809 m_singularValues.coeffRef(i) = a;
810 if(computeU() && (a!=
RealScalar(0))) m_matrixU.col(i) *= m_workMatrix.coeff(i,i)/a;
815 m_nonzeroSingularValues = m_diagSize;
816 for(
Index i = 0; i < m_diagSize; i++)
819 RealScalar maxRemainingSingularValue = m_singularValues.tail(m_diagSize-i).maxCoeff(&pos);
820 if(maxRemainingSingularValue ==
RealScalar(0))
822 m_nonzeroSingularValues = i;
828 std::swap(m_singularValues.coeffRef(i), m_singularValues.coeffRef(pos));
829 if(computeU()) m_matrixU.col(pos).swap(m_matrixU.col(i));
830 if(computeV()) m_matrixV.col(pos).swap(m_matrixV.col(i));
834 m_isInitialized =
true;
839 template<
typename _MatrixType,
int QRPreconditioner,
typename Rhs>
846 template<typename Dest>
void evalTo(Dest& dst)
const 854 Index nonzeroSingVals = dec().nonzeroSingularValues();
856 tmp.noalias() = dec().matrixU().leftCols(nonzeroSingVals).adjoint() * rhs();
857 tmp = dec().singularValues().head(nonzeroSingVals).asDiagonal().inverse() * tmp;
858 dst = dec().matrixV().leftCols(nonzeroSingVals) * tmp;
870 template<
typename Derived>
879 #endif // EIGEN_JACOBISVD_H
Matrix< Scalar, DiagSizeAtCompileTime, DiagSizeAtCompileTime, MatrixOptions, MaxDiagSizeAtCompileTime, MaxDiagSizeAtCompileTime > WorkMatrixType
Matrix< Scalar, ColsAtCompileTime, RowsAtCompileTime, Options, MaxColsAtCompileTime, MaxRowsAtCompileTime > TransposeTypeWithSameStorageOrder
bool run(JacobiSVD< MatrixType, QRPreconditioner > &, const MatrixType &)
MatrixType::Scalar Scalar
internal::plain_col_type< MatrixType >::type ColType
Matrix< Scalar, ColsAtCompileTime, RowsAtCompileTime, Options, MaxColsAtCompileTime, MaxRowsAtCompileTime > TransposeTypeWithSameStorageOrder
void allocate(const JacobiSVD< MatrixType, QRPreconditioner > &)
Index nonzeroSingularValues() const
JacobiSVD< _MatrixType, QRPreconditioner > JacobiSVDType
unsigned int m_computationOptions
RealReturnType real() const
static void run(typename SVD::WorkMatrixType &work_matrix, SVD &svd, Index p, Index q)
internal::plain_col_type< MatrixType >::type m_workspace
void allocate(const JacobiSVD< MatrixType, HouseholderQRPreconditioner > &svd)
Rotation given by a cosine-sine pair.
MatrixType::Scalar Scalar
void allocate(const JacobiSVD< MatrixType, FullPivHouseholderQRPreconditioner > &svd)
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
TransposeTypeWithSameStorageOrder m_adjoint
#define EIGEN_SIZE_MIN_PREFER_FIXED(a, b)
internal::plain_diag_type< MatrixType, RealScalar >::type SingularValuesType
MatrixType::Scalar Scalar
WorkspaceType m_workspace
const MatrixUType & matrixU() const
bool run(JacobiSVD< MatrixType, ColPivHouseholderQRPreconditioner > &svd, const MatrixType &matrix)
void real_2x2_jacobi_svd(const MatrixType &matrix, Index p, Index q, JacobiRotation< RealScalar > *j_left, JacobiRotation< RealScalar > *j_right)
TransposeTypeWithSameStorageOrder m_adjoint
FullPivHouseholderQR< MatrixType > QRType
MatrixType::Scalar Scalar
JacobiSVD< PlainObject > jacobiSvd(unsigned int computationOptions=0) const
static void run(typename SVD::WorkMatrixType &, SVD &, Index, Index)
JacobiSVD()
Default Constructor.
const internal::solve_retval< JacobiSVD, Rhs > solve(const MatrixBase< Rhs > &b) const
bool makeJacobi(const MatrixBase< Derived > &, typename Derived::Index p, typename Derived::Index q)
internal::qr_preconditioner_impl< MatrixType, QRPreconditioner, internal::PreconditionIfMoreColsThanRows > m_qr_precond_morecols
HouseholderQR< MatrixType > QRType
ColPivHouseholderQR< TransposeTypeWithSameStorageOrder > QRType
Matrix< Scalar, ColsAtCompileTime, RowsAtCompileTime, Options, MaxColsAtCompileTime, MaxRowsAtCompileTime > TransposeTypeWithSameStorageOrder
Matrix< Scalar, RowsAtCompileTime, RowsAtCompileTime, MatrixOptions, MaxRowsAtCompileTime, MaxRowsAtCompileTime > MatrixUType
internal::plain_row_type< MatrixType >::type m_workspace
FullPivHouseholderQR< TransposeTypeWithSameStorageOrder > QRType
SingularValuesType m_singularValues
TransposeTypeWithSameStorageOrder m_adjoint
const MatrixVType & matrixV() const
bool run(JacobiSVD< MatrixType, HouseholderQRPreconditioner > &svd, const MatrixType &matrix)
bool run(JacobiSVD< MatrixType, FullPivHouseholderQRPreconditioner > &svd, const MatrixType &matrix)
MatrixType::RealScalar RealScalar
Matrix< Scalar, 1, RowsAtCompileTime, RowMajor, 1, MaxRowsAtCompileTime > WorkspaceType
bool run(JacobiSVD< MatrixType, HouseholderQRPreconditioner > &svd, const MatrixType &matrix)
void allocate(const JacobiSVD< MatrixType, ColPivHouseholderQRPreconditioner > &svd)
void allocate(const JacobiSVD< MatrixType, HouseholderQRPreconditioner > &svd)
ColPivHouseholderQR< MatrixType > QRType
HouseholderQR< TransposeTypeWithSameStorageOrder > QRType
EIGEN_STRONG_INLINE const CwiseUnaryOp< internal::scalar_abs_op< Scalar >, const Derived > abs() const
#define EIGEN_SIZE_MIN_PREFER_DYNAMIC(a, b)
JacobiSVD & compute(const MatrixType &matrix)
Method performing the decomposition of given matrix using current options.
Matrix< Scalar, ColsAtCompileTime, ColsAtCompileTime, MatrixOptions, MaxColsAtCompileTime, MaxColsAtCompileTime > MatrixVType
void allocate(Index rows, Index cols, unsigned int computationOptions)
Two-sided Jacobi SVD decomposition of a rectangular matrix.
JacobiRotation transpose() const
JacobiSVD< MatrixType, QRPreconditioner > SVD
MatrixType::Scalar Scalar
internal::plain_row_type< MatrixType >::type RowType
NumTraits< typename MatrixType::Scalar >::Real RealScalar
EIGEN_STRONG_INLINE const Scalar & coeff(Index rowId, Index colId) const
JacobiSVD< MatrixType, QRPreconditioner > SVD
JacobiSVD & compute(const MatrixType &matrix, unsigned int computationOptions)
Method performing the decomposition of given matrix using custom options.
bool run(JacobiSVD< MatrixType, FullPivHouseholderQRPreconditioner > &svd, const MatrixType &matrix)
internal::plain_col_type< MatrixType >::type m_workspace
#define EIGEN_MAKE_SOLVE_HELPERS(DecompositionType, Rhs)
#define EIGEN_IMPLIES(a, b)
internal::qr_preconditioner_impl< MatrixType, QRPreconditioner, internal::PreconditionIfMoreRowsThanCols > m_qr_precond_morerows
const SingularValuesType & singularValues() const
internal::plain_row_type< MatrixType >::type m_workspace
EIGEN_STRONG_INLINE const CwiseUnaryOp< internal::scalar_abs2_op< Scalar >, const Derived > abs2() const
JacobiSVD(Index rows, Index cols, unsigned int computationOptions=0)
Default Constructor with memory preallocation.
Base class for all dense matrices, vectors, and expressions.
bool run(JacobiSVD< MatrixType, ColPivHouseholderQRPreconditioner > &svd, const MatrixType &matrix)
WorkMatrixType m_workMatrix
JacobiSVD(const MatrixType &matrix, unsigned int computationOptions=0)
Constructor performing the decomposition of given matrix.
JacobiRotation adjoint() const
internal::plain_row_type< MatrixType >::type m_workspace
void allocate(const JacobiSVD< MatrixType, ColPivHouseholderQRPreconditioner > &svd)
void allocate(const JacobiSVD< MatrixType, FullPivHouseholderQRPreconditioner > &svd)
MatrixType::Scalar Scalar
const CwiseUnaryOp< internal::scalar_sqrt_op< Scalar >, const Derived > sqrt() const