10 #ifndef EIGEN_JACOBISVD_H
11 #define EIGEN_JACOBISVD_H
18 template<
typename MatrixType,
int QRPreconditioner,
19 bool IsComplex = NumTraits<typename MatrixType::Scalar>::IsComplex>
20 struct svd_precondition_2x2_block_to_be_real {};
29 enum { PreconditionIfMoreColsThanRows, PreconditionIfMoreRowsThanCols };
31 template<
typename MatrixType,
int QRPreconditioner,
int Case>
32 struct qr_preconditioner_should_do_anything
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,
40 ret = !( (QRPreconditioner == NoQRPreconditioner) ||
41 (Case == PreconditionIfMoreColsThanRows &&
bool(a)) ||
42 (Case == PreconditionIfMoreRowsThanCols &&
bool(b)) )
46 template<
typename MatrixType,
int QRPreconditioner,
int Case,
47 bool DoAnything = qr_preconditioner_should_do_anything<MatrixType, QRPreconditioner, Case>::ret
48 >
struct qr_preconditioner_impl {};
50 template<
typename MatrixType,
int QRPreconditioner,
int Case>
51 class qr_preconditioner_impl<MatrixType, QRPreconditioner, Case, false>
54 typedef typename MatrixType::Index Index;
55 void allocate(
const JacobiSVD<MatrixType, QRPreconditioner>&) {}
56 bool run(JacobiSVD<MatrixType, QRPreconditioner>&,
const MatrixType&)
64 template<
typename MatrixType>
65 class qr_preconditioner_impl<MatrixType, FullPivHouseholderQRPreconditioner, PreconditionIfMoreRowsThanCols, true>
68 typedef typename MatrixType::Index Index;
69 typedef typename MatrixType::Scalar Scalar;
72 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
73 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime
75 typedef Matrix<Scalar, 1, RowsAtCompileTime, RowMajor, 1, MaxRowsAtCompileTime> WorkspaceType;
77 void allocate(
const JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner>& svd)
79 if (svd.rows() != m_qr.rows() || svd.cols() != m_qr.cols())
82 ::new (&m_qr) QRType(svd.rows(), svd.cols());
84 if (svd.m_computeFullU) m_workspace.resize(svd.rows());
87 bool run(JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
89 if(matrix.rows() > matrix.cols())
92 svd.m_workMatrix = m_qr.matrixQR().block(0,0,matrix.cols(),matrix.cols()).
template triangularView<Upper>();
93 if(svd.m_computeFullU) m_qr.matrixQ().evalTo(svd.m_matrixU, m_workspace);
94 if(svd.computeV()) svd.m_matrixV = m_qr.colsPermutation();
100 typedef FullPivHouseholderQR<MatrixType> QRType;
102 WorkspaceType m_workspace;
105 template<
typename MatrixType>
106 class qr_preconditioner_impl<MatrixType, FullPivHouseholderQRPreconditioner, PreconditionIfMoreColsThanRows, true>
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
119 typedef Matrix<Scalar, ColsAtCompileTime, RowsAtCompileTime, Options, MaxColsAtCompileTime, MaxRowsAtCompileTime>
120 TransposeTypeWithSameStorageOrder;
122 void allocate(
const JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner>& svd)
124 if (svd.cols() != m_qr.rows() || svd.rows() != m_qr.cols())
127 ::new (&m_qr) QRType(svd.cols(), svd.rows());
129 m_adjoint.resize(svd.cols(), svd.rows());
130 if (svd.m_computeFullV) m_workspace.resize(svd.cols());
133 bool run(JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
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();
140 if(svd.m_computeFullV) m_qr.matrixQ().evalTo(svd.m_matrixV, m_workspace);
141 if(svd.computeU()) svd.m_matrixU = m_qr.colsPermutation();
147 typedef FullPivHouseholderQR<TransposeTypeWithSameStorageOrder> QRType;
149 TransposeTypeWithSameStorageOrder m_adjoint;
150 typename internal::plain_row_type<MatrixType>::type m_workspace;
155 template<
typename MatrixType>
156 class qr_preconditioner_impl<MatrixType, ColPivHouseholderQRPreconditioner, PreconditionIfMoreRowsThanCols, true>
159 typedef typename MatrixType::Index Index;
161 void allocate(
const JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner>& svd)
163 if (svd.rows() != m_qr.rows() || svd.cols() != m_qr.cols())
166 ::new (&m_qr) QRType(svd.rows(), svd.cols());
168 if (svd.m_computeFullU) m_workspace.resize(svd.rows());
169 else if (svd.m_computeThinU) m_workspace.resize(svd.cols());
172 bool run(JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
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>();
178 if(svd.m_computeFullU) m_qr.householderQ().evalTo(svd.m_matrixU, m_workspace);
179 else if(svd.m_computeThinU)
181 svd.m_matrixU.setIdentity(matrix.rows(), matrix.cols());
182 m_qr.householderQ().applyThisOnTheLeft(svd.m_matrixU, m_workspace);
184 if(svd.computeV()) svd.m_matrixV = m_qr.colsPermutation();
191 typedef ColPivHouseholderQR<MatrixType> QRType;
193 typename internal::plain_col_type<MatrixType>::type m_workspace;
196 template<
typename MatrixType>
197 class qr_preconditioner_impl<MatrixType, ColPivHouseholderQRPreconditioner, PreconditionIfMoreColsThanRows, true>
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
211 typedef Matrix<Scalar, ColsAtCompileTime, RowsAtCompileTime, Options, MaxColsAtCompileTime, MaxRowsAtCompileTime>
212 TransposeTypeWithSameStorageOrder;
214 void allocate(
const JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner>& svd)
216 if (svd.cols() != m_qr.rows() || svd.rows() != m_qr.cols())
219 ::new (&m_qr) QRType(svd.cols(), svd.rows());
221 if (svd.m_computeFullV) m_workspace.resize(svd.cols());
222 else if (svd.m_computeThinV) m_workspace.resize(svd.rows());
223 m_adjoint.resize(svd.cols(), svd.rows());
226 bool run(JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
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();
234 if(svd.m_computeFullV) m_qr.householderQ().evalTo(svd.m_matrixV, m_workspace);
235 else if(svd.m_computeThinV)
237 svd.m_matrixV.setIdentity(matrix.cols(), matrix.rows());
238 m_qr.householderQ().applyThisOnTheLeft(svd.m_matrixV, m_workspace);
240 if(svd.computeU()) svd.m_matrixU = m_qr.colsPermutation();
247 typedef ColPivHouseholderQR<TransposeTypeWithSameStorageOrder> QRType;
249 TransposeTypeWithSameStorageOrder m_adjoint;
250 typename internal::plain_row_type<MatrixType>::type m_workspace;
255 template<
typename MatrixType>
256 class qr_preconditioner_impl<MatrixType, HouseholderQRPreconditioner, PreconditionIfMoreRowsThanCols, true>
259 typedef typename MatrixType::Index Index;
261 void allocate(
const JacobiSVD<MatrixType, HouseholderQRPreconditioner>& svd)
263 if (svd.rows() != m_qr.rows() || svd.cols() != m_qr.cols())
266 ::new (&m_qr) QRType(svd.rows(), svd.cols());
268 if (svd.m_computeFullU) m_workspace.resize(svd.rows());
269 else if (svd.m_computeThinU) m_workspace.resize(svd.cols());
272 bool run(JacobiSVD<MatrixType, HouseholderQRPreconditioner>& svd, const MatrixType& matrix)
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>();
278 if(svd.m_computeFullU) m_qr.householderQ().evalTo(svd.m_matrixU, m_workspace);
279 else if(svd.m_computeThinU)
281 svd.m_matrixU.setIdentity(matrix.rows(), matrix.cols());
282 m_qr.householderQ().applyThisOnTheLeft(svd.m_matrixU, m_workspace);
284 if(svd.computeV()) svd.m_matrixV.setIdentity(matrix.cols(), matrix.cols());
290 typedef HouseholderQR<MatrixType> QRType;
292 typename internal::plain_col_type<MatrixType>::type m_workspace;
295 template<
typename MatrixType>
296 class qr_preconditioner_impl<MatrixType, HouseholderQRPreconditioner, PreconditionIfMoreColsThanRows, true>
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
310 typedef Matrix<Scalar, ColsAtCompileTime, RowsAtCompileTime, Options, MaxColsAtCompileTime, MaxRowsAtCompileTime>
311 TransposeTypeWithSameStorageOrder;
313 void allocate(
const JacobiSVD<MatrixType, HouseholderQRPreconditioner>& svd)
315 if (svd.cols() != m_qr.rows() || svd.rows() != m_qr.cols())
318 ::new (&m_qr) QRType(svd.cols(), svd.rows());
320 if (svd.m_computeFullV) m_workspace.resize(svd.cols());
321 else if (svd.m_computeThinV) m_workspace.resize(svd.rows());
322 m_adjoint.resize(svd.cols(), svd.rows());
325 bool run(JacobiSVD<MatrixType, HouseholderQRPreconditioner>& svd, const MatrixType& matrix)
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();
333 if(svd.m_computeFullV) m_qr.householderQ().evalTo(svd.m_matrixV, m_workspace);
334 else if(svd.m_computeThinV)
336 svd.m_matrixV.setIdentity(matrix.cols(), matrix.rows());
337 m_qr.householderQ().applyThisOnTheLeft(svd.m_matrixV, m_workspace);
339 if(svd.computeU()) svd.m_matrixU.setIdentity(matrix.rows(), matrix.rows());
346 typedef HouseholderQR<TransposeTypeWithSameStorageOrder> QRType;
348 TransposeTypeWithSameStorageOrder m_adjoint;
349 typename internal::plain_row_type<MatrixType>::type m_workspace;
357 template<
typename MatrixType,
int QRPreconditioner>
358 struct svd_precondition_2x2_block_to_be_real<MatrixType, QRPreconditioner, false>
360 typedef JacobiSVD<MatrixType, QRPreconditioner> SVD;
361 typedef typename SVD::Index Index;
362 static void run(
typename SVD::WorkMatrixType&, SVD&, Index, Index) {}
365 template<
typename MatrixType,
int QRPreconditioner>
366 struct svd_precondition_2x2_block_to_be_real<MatrixType, QRPreconditioner, true>
368 typedef JacobiSVD<MatrixType, QRPreconditioner> SVD;
369 typedef typename MatrixType::Scalar Scalar;
370 typedef typename MatrixType::RealScalar RealScalar;
371 typedef typename SVD::Index Index;
372 static void run(
typename SVD::WorkMatrixType& work_matrix, SVD& svd, Index p, Index q)
376 JacobiRotation<Scalar> rot;
377 RealScalar n = sqrt(numext::abs2(work_matrix.coeff(p,p)) + numext::abs2(work_matrix.coeff(q,p)));
380 z = abs(work_matrix.coeff(p,q)) / work_matrix.coeff(p,q);
381 work_matrix.row(p) *= z;
382 if(svd.computeU()) svd.m_matrixU.col(p) *= conj(z);
383 z = abs(work_matrix.coeff(q,q)) / work_matrix.coeff(q,q);
384 work_matrix.row(q) *= z;
385 if(svd.computeU()) svd.m_matrixU.col(q) *= conj(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);
392 if(svd.computeU()) svd.m_matrixU.applyOnTheRight(p,q,rot.adjoint());
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;
397 if(svd.computeV()) svd.m_matrixV.col(q) *= z;
399 if(work_matrix.coeff(q,q) != Scalar(0))
401 z = abs(work_matrix.coeff(q,q)) / work_matrix.coeff(q,q);
402 work_matrix.row(q) *= z;
403 if(svd.computeU()) svd.m_matrixU.col(q) *= conj(z);
409 template<
typename MatrixType,
typename RealScalar,
typename Index>
410 void real_2x2_jacobi_svd(
const MatrixType& matrix, Index p, Index q,
411 JacobiRotation<RealScalar> *j_left,
412 JacobiRotation<RealScalar> *j_right)
415 Matrix<RealScalar,2,2> m;
416 m << numext::real(matrix.coeff(p,p)), numext::real(matrix.coeff(p,q)),
417 numext::real(matrix.coeff(q,p)), numext::real(matrix.coeff(q,q));
418 JacobiRotation<RealScalar> rot1;
419 RealScalar t = m.coeff(0,0) + m.coeff(1,1);
420 RealScalar d = m.coeff(1,0) - m.coeff(0,1);
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;
429 rot1.c() = RealScalar(1) / sqrt(RealScalar(1) + numext::abs2(u));
430 rot1.s() = rot1.c() * u;
432 m.applyOnTheLeft(0,1,rot1);
433 j_right->makeJacobi(m,0,1);
434 *j_left = rot1 * j_right->transpose();
492 template<
typename _MatrixType,
int QRPreconditioner>
497 typedef _MatrixType MatrixType;
498 typedef typename MatrixType::Scalar Scalar;
499 typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
500 typedef typename MatrixType::Index Index;
502 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
503 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
504 DiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_DYNAMIC(RowsAtCompileTime,ColsAtCompileTime),
505 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
506 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
507 MaxDiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(MaxRowsAtCompileTime,MaxColsAtCompileTime),
508 MatrixOptions = MatrixType::Options
511 typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime,
512 MatrixOptions, MaxRowsAtCompileTime, MaxRowsAtCompileTime>
514 typedef Matrix<Scalar, ColsAtCompileTime, ColsAtCompileTime,
515 MatrixOptions, MaxColsAtCompileTime, MaxColsAtCompileTime>
517 typedef typename internal::plain_diag_type<MatrixType, RealScalar>::type SingularValuesType;
518 typedef typename internal::plain_row_type<MatrixType>::type RowType;
519 typedef typename internal::plain_col_type<MatrixType>::type ColType;
520 typedef Matrix<Scalar, DiagSizeAtCompileTime, DiagSizeAtCompileTime,
521 MatrixOptions, MaxDiagSizeAtCompileTime, MaxDiagSizeAtCompileTime>
540 JacobiSVD(Index rows, Index cols,
unsigned int computationOptions = 0)
543 allocate(rows, cols, computationOptions);
556 JacobiSVD(
const MatrixType& matrix,
unsigned int computationOptions = 0)
559 compute(matrix, computationOptions);
582 return compute(matrix, this->m_computationOptions);
594 template<
typename Rhs>
595 inline const internal::solve_retval<JacobiSVD, Rhs>
596 solve(
const MatrixBase<Rhs>& b)
const
598 eigen_assert(this->m_isInitialized &&
"JacobiSVD is not initialized.");
600 return internal::solve_retval<JacobiSVD, Rhs>(*
this, b.derived());
606 void allocate(Index rows, Index cols,
unsigned int computationOptions);
609 WorkMatrixType m_workMatrix;
611 template<
typename __MatrixType,
int _QRPreconditioner,
bool _IsComplex>
612 friend struct internal::svd_precondition_2x2_block_to_be_real;
613 template<
typename __MatrixType,
int _QRPreconditioner,
int _Case,
bool _DoAnything>
614 friend struct internal::qr_preconditioner_impl;
616 internal::qr_preconditioner_impl<MatrixType, QRPreconditioner, internal::PreconditionIfMoreColsThanRows> m_qr_precond_morecols;
617 internal::qr_preconditioner_impl<MatrixType, QRPreconditioner, internal::PreconditionIfMoreRowsThanCols> m_qr_precond_morerows;
620 template<
typename MatrixType,
int QRPreconditioner>
621 void JacobiSVD<MatrixType, QRPreconditioner>::allocate(Index rows, Index cols,
unsigned int computationOptions)
623 if (SVDBase<MatrixType>::allocate(rows, cols, computationOptions))
return;
625 if (QRPreconditioner == FullPivHouseholderQRPreconditioner)
627 eigen_assert(!(this->m_computeThinU || this->m_computeThinV) &&
628 "JacobiSVD: can't compute thin U or thin V with the FullPivHouseholderQR preconditioner. "
629 "Use the ColPivHouseholderQR preconditioner instead.");
632 m_workMatrix.resize(this->m_diagSize, this->m_diagSize);
634 if(this->m_cols>this->m_rows) m_qr_precond_morecols.allocate(*
this);
635 if(this->m_rows>this->m_cols) m_qr_precond_morerows.allocate(*
this);
638 template<
typename MatrixType,
int QRPreconditioner>
643 allocate(matrix.rows(), matrix.cols(), computationOptions);
647 const RealScalar precision = RealScalar(2) * NumTraits<Scalar>::epsilon();
650 const RealScalar considerAsZero = RealScalar(2) * std::numeric_limits<RealScalar>::denorm_min();
654 if(!m_qr_precond_morecols.run(*
this, matrix) && !m_qr_precond_morerows.run(*
this, matrix))
656 m_workMatrix = matrix.block(0,0,this->m_diagSize,this->m_diagSize);
657 if(this->m_computeFullU) this->m_matrixU.setIdentity(this->m_rows,this->m_rows);
658 if(this->m_computeThinU) this->m_matrixU.setIdentity(this->m_rows,this->m_diagSize);
659 if(this->m_computeFullV) this->m_matrixV.setIdentity(this->m_cols,this->m_cols);
660 if(this->m_computeThinV) this->m_matrixV.setIdentity(this->m_cols, this->m_diagSize);
665 bool finished =
false;
672 for(Index p = 1; p < this->m_diagSize; ++p)
674 for(Index q = 0; q < p; ++q)
680 RealScalar threshold = (max)(considerAsZero, precision * (max)(abs(m_workMatrix.coeff(p,p)),
681 abs(m_workMatrix.coeff(q,q))));
682 if((max)(abs(m_workMatrix.coeff(p,q)),abs(m_workMatrix.coeff(q,p))) > threshold)
687 internal::svd_precondition_2x2_block_to_be_real<MatrixType, QRPreconditioner>::run(m_workMatrix, *
this, p, q);
688 JacobiRotation<RealScalar> j_left, j_right;
689 internal::real_2x2_jacobi_svd(m_workMatrix, p, q, &j_left, &j_right);
692 m_workMatrix.applyOnTheLeft(p,q,j_left);
695 m_workMatrix.applyOnTheRight(p,q,j_right);
704 for(Index i = 0; i < this->m_diagSize; ++i)
706 RealScalar a = abs(m_workMatrix.coeff(i,i));
707 this->m_singularValues.coeffRef(i) = a;
713 this->m_nonzeroSingularValues = this->m_diagSize;
714 for(Index i = 0; i < this->m_diagSize; i++)
717 RealScalar maxRemainingSingularValue = this->m_singularValues.tail(this->m_diagSize-i).maxCoeff(&pos);
718 if(maxRemainingSingularValue == RealScalar(0))
720 this->m_nonzeroSingularValues = i;
726 std::swap(this->m_singularValues.coeffRef(i), this->m_singularValues.coeffRef(pos));
732 this->m_isInitialized =
true;
737 template<
typename _MatrixType,
int QRPreconditioner,
typename Rhs>
738 struct solve_retval<
JacobiSVD<_MatrixType, QRPreconditioner>, Rhs>
739 : solve_retval_base<JacobiSVD<_MatrixType, QRPreconditioner>, Rhs>
742 EIGEN_MAKE_SOLVE_HELPERS(JacobiSVDType,Rhs)
744 template<typename Dest>
void evalTo(Dest& dst)
const
746 eigen_assert(rhs().rows() == dec().rows());
751 Index diagSize = (std::min)(dec().rows(), dec().cols());
752 typename JacobiSVDType::SingularValuesType invertedSingVals(diagSize);
754 Index nonzeroSingVals = dec().nonzeroSingularValues();
755 invertedSingVals.head(nonzeroSingVals) = dec().singularValues().head(nonzeroSingVals).array().inverse();
756 invertedSingVals.tail(diagSize - nonzeroSingVals).setZero();
758 dst = dec().matrixV().leftCols(diagSize)
759 * invertedSingVals.asDiagonal()
760 * dec().matrixU().leftCols(diagSize).adjoint()
773 template<
typename Derived>
774 JacobiSVD<typename MatrixBase<Derived>::PlainObject>
775 MatrixBase<Derived>::jacobiSvd(
unsigned int computationOptions)
const
777 return JacobiSVD<PlainObject>(*
this, computationOptions);
782 #endif // EIGEN_JACOBISVD_H
const internal::solve_retval< JacobiSVD, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: JacobiSVD.h:596
SVDBase< MatrixType > & compute(const MatrixType &matrix)
Method performing the decomposition of given matrix using current options.
Definition: JacobiSVD.h:580
SVDBase< MatrixType > & compute(const MatrixType &matrix, unsigned int computationOptions)
Method performing the decomposition of given matrix using custom options.
Definition: JacobiSVD.h:640
JacobiSVD()
Default Constructor.
Definition: JacobiSVD.h:529
Mother class of SVD classes algorithms.
Definition: SVDBase.h:46
JacobiSVD(Index rows, Index cols, unsigned int computationOptions=0)
Default Constructor with memory preallocation.
Definition: JacobiSVD.h:540
Two-sided Jacobi SVD decomposition of a rectangular matrix.
Definition: JacobiSVD.h:493
JacobiSVD(const MatrixType &matrix, unsigned int computationOptions=0)
Constructor performing the decomposition of given matrix.
Definition: JacobiSVD.h:556