1
// This file is part of Eigen, a lightweight C++ template library
4
// Copyright (C) 2009-2010 Benoit Jacob <jacob.benoit.1@gmail.com>
6
// Eigen is free software; you can redistribute it and/or
7
// modify it under the terms of the GNU Lesser General Public
8
// License as published by the Free Software Foundation; either
9
// version 3 of the License, or (at your option) any later version.
11
// Alternatively, you can redistribute it and/or
12
// modify it under the terms of the GNU General Public License as
13
// published by the Free Software Foundation; either version 2 of
14
// the License, or (at your option) any later version.
16
// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
17
// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
18
// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
19
// GNU General Public License for more details.
21
// You should have received a copy of the GNU Lesser General Public
22
// License and a copy of the GNU General Public License along with
23
// Eigen. If not, see <http://www.gnu.org/licenses/>.
25
#ifndef EIGEN_JACOBISVD_H
26
#define EIGEN_JACOBISVD_H
29
// forward declaration (needed by ICC)
30
// the empty body is required by MSVC
31
template<typename MatrixType, int QRPreconditioner,
32
bool IsComplex = NumTraits<typename MatrixType::Scalar>::IsComplex>
33
struct svd_precondition_2x2_block_to_be_real {};
35
/*** QR preconditioners (R-SVD)
37
*** Their role is to reduce the problem of computing the SVD to the case of a square matrix.
38
*** This approach, known as R-SVD, is an optimization for rectangular-enough matrices, and is a requirement for
39
*** JacobiSVD which by itself is only able to work on square matrices.
42
enum { PreconditionIfMoreColsThanRows, PreconditionIfMoreRowsThanCols };
44
template<typename MatrixType, int QRPreconditioner, int Case>
45
struct qr_preconditioner_should_do_anything
47
enum { a = MatrixType::RowsAtCompileTime != Dynamic &&
48
MatrixType::ColsAtCompileTime != Dynamic &&
49
MatrixType::ColsAtCompileTime <= MatrixType::RowsAtCompileTime,
50
b = MatrixType::RowsAtCompileTime != Dynamic &&
51
MatrixType::ColsAtCompileTime != Dynamic &&
52
MatrixType::RowsAtCompileTime <= MatrixType::ColsAtCompileTime,
53
ret = !( (QRPreconditioner == NoQRPreconditioner) ||
54
(Case == PreconditionIfMoreColsThanRows && bool(a)) ||
55
(Case == PreconditionIfMoreRowsThanCols && bool(b)) )
59
template<typename MatrixType, int QRPreconditioner, int Case,
60
bool DoAnything = qr_preconditioner_should_do_anything<MatrixType, QRPreconditioner, Case>::ret
61
> struct qr_preconditioner_impl {};
63
template<typename MatrixType, int QRPreconditioner, int Case>
64
struct qr_preconditioner_impl<MatrixType, QRPreconditioner, Case, false>
66
static bool run(JacobiSVD<MatrixType, QRPreconditioner>&, const MatrixType&)
72
/*** preconditioner using FullPivHouseholderQR ***/
74
template<typename MatrixType>
75
struct qr_preconditioner_impl<MatrixType, FullPivHouseholderQRPreconditioner, PreconditionIfMoreRowsThanCols, true>
77
static bool run(JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
79
if(matrix.rows() > matrix.cols())
81
FullPivHouseholderQR<MatrixType> qr(matrix);
82
svd.m_workMatrix = qr.matrixQR().block(0,0,matrix.cols(),matrix.cols()).template triangularView<Upper>();
83
if(svd.m_computeFullU) svd.m_matrixU = qr.matrixQ();
84
if(svd.computeV()) svd.m_matrixV = qr.colsPermutation();
91
template<typename MatrixType>
92
struct qr_preconditioner_impl<MatrixType, FullPivHouseholderQRPreconditioner, PreconditionIfMoreColsThanRows, true>
94
static bool run(JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
96
if(matrix.cols() > matrix.rows())
98
typedef Matrix<typename MatrixType::Scalar, MatrixType::ColsAtCompileTime, MatrixType::RowsAtCompileTime,
99
MatrixType::Options, MatrixType::MaxColsAtCompileTime, MatrixType::MaxRowsAtCompileTime>
100
TransposeTypeWithSameStorageOrder;
101
FullPivHouseholderQR<TransposeTypeWithSameStorageOrder> qr(matrix.adjoint());
102
svd.m_workMatrix = qr.matrixQR().block(0,0,matrix.rows(),matrix.rows()).template triangularView<Upper>().adjoint();
103
if(svd.m_computeFullV) svd.m_matrixV = qr.matrixQ();
104
if(svd.computeU()) svd.m_matrixU = qr.colsPermutation();
111
/*** preconditioner using ColPivHouseholderQR ***/
113
template<typename MatrixType>
114
struct qr_preconditioner_impl<MatrixType, ColPivHouseholderQRPreconditioner, PreconditionIfMoreRowsThanCols, true>
116
static bool run(JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
118
if(matrix.rows() > matrix.cols())
120
ColPivHouseholderQR<MatrixType> qr(matrix);
121
svd.m_workMatrix = qr.matrixQR().block(0,0,matrix.cols(),matrix.cols()).template triangularView<Upper>();
122
if(svd.m_computeFullU) svd.m_matrixU = qr.householderQ();
123
else if(svd.m_computeThinU) {
124
svd.m_matrixU.setIdentity(matrix.rows(), matrix.cols());
125
qr.householderQ().applyThisOnTheLeft(svd.m_matrixU);
127
if(svd.computeV()) svd.m_matrixV = qr.colsPermutation();
134
template<typename MatrixType>
135
struct qr_preconditioner_impl<MatrixType, ColPivHouseholderQRPreconditioner, PreconditionIfMoreColsThanRows, true>
137
static bool run(JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
139
if(matrix.cols() > matrix.rows())
141
typedef Matrix<typename MatrixType::Scalar, MatrixType::ColsAtCompileTime, MatrixType::RowsAtCompileTime,
142
MatrixType::Options, MatrixType::MaxColsAtCompileTime, MatrixType::MaxRowsAtCompileTime>
143
TransposeTypeWithSameStorageOrder;
144
ColPivHouseholderQR<TransposeTypeWithSameStorageOrder> qr(matrix.adjoint());
145
svd.m_workMatrix = qr.matrixQR().block(0,0,matrix.rows(),matrix.rows()).template triangularView<Upper>().adjoint();
146
if(svd.m_computeFullV) svd.m_matrixV = qr.householderQ();
147
else if(svd.m_computeThinV) {
148
svd.m_matrixV.setIdentity(matrix.cols(), matrix.rows());
149
qr.householderQ().applyThisOnTheLeft(svd.m_matrixV);
151
if(svd.computeU()) svd.m_matrixU = qr.colsPermutation();
158
/*** preconditioner using HouseholderQR ***/
160
template<typename MatrixType>
161
struct qr_preconditioner_impl<MatrixType, HouseholderQRPreconditioner, PreconditionIfMoreRowsThanCols, true>
163
static bool run(JacobiSVD<MatrixType, HouseholderQRPreconditioner>& svd, const MatrixType& matrix)
165
if(matrix.rows() > matrix.cols())
167
HouseholderQR<MatrixType> qr(matrix);
168
svd.m_workMatrix = qr.matrixQR().block(0,0,matrix.cols(),matrix.cols()).template triangularView<Upper>();
169
if(svd.m_computeFullU) svd.m_matrixU = qr.householderQ();
170
else if(svd.m_computeThinU) {
171
svd.m_matrixU.setIdentity(matrix.rows(), matrix.cols());
172
qr.householderQ().applyThisOnTheLeft(svd.m_matrixU);
174
if(svd.computeV()) svd.m_matrixV.setIdentity(matrix.cols(), matrix.cols());
181
template<typename MatrixType>
182
struct qr_preconditioner_impl<MatrixType, HouseholderQRPreconditioner, PreconditionIfMoreColsThanRows, true>
184
static bool run(JacobiSVD<MatrixType, HouseholderQRPreconditioner>& svd, const MatrixType& matrix)
186
if(matrix.cols() > matrix.rows())
188
typedef Matrix<typename MatrixType::Scalar, MatrixType::ColsAtCompileTime, MatrixType::RowsAtCompileTime,
189
MatrixType::Options, MatrixType::MaxColsAtCompileTime, MatrixType::MaxRowsAtCompileTime>
190
TransposeTypeWithSameStorageOrder;
191
HouseholderQR<TransposeTypeWithSameStorageOrder> qr(matrix.adjoint());
192
svd.m_workMatrix = qr.matrixQR().block(0,0,matrix.rows(),matrix.rows()).template triangularView<Upper>().adjoint();
193
if(svd.m_computeFullV) svd.m_matrixV = qr.householderQ();
194
else if(svd.m_computeThinV) {
195
svd.m_matrixV.setIdentity(matrix.cols(), matrix.rows());
196
qr.householderQ().applyThisOnTheLeft(svd.m_matrixV);
198
if(svd.computeU()) svd.m_matrixU.setIdentity(matrix.rows(), matrix.rows());
205
/*** 2x2 SVD implementation
207
*** JacobiSVD consists in performing a series of 2x2 SVD subproblems
210
template<typename MatrixType, int QRPreconditioner>
211
struct svd_precondition_2x2_block_to_be_real<MatrixType, QRPreconditioner, false>
213
typedef JacobiSVD<MatrixType, QRPreconditioner> SVD;
214
typedef typename SVD::Index Index;
215
static void run(typename SVD::WorkMatrixType&, SVD&, Index, Index) {}
218
template<typename MatrixType, int QRPreconditioner>
219
struct svd_precondition_2x2_block_to_be_real<MatrixType, QRPreconditioner, true>
221
typedef JacobiSVD<MatrixType, QRPreconditioner> SVD;
222
typedef typename MatrixType::Scalar Scalar;
223
typedef typename MatrixType::RealScalar RealScalar;
224
typedef typename SVD::Index Index;
225
static void run(typename SVD::WorkMatrixType& work_matrix, SVD& svd, Index p, Index q)
228
JacobiRotation<Scalar> rot;
229
RealScalar n = sqrt(abs2(work_matrix.coeff(p,p)) + abs2(work_matrix.coeff(q,p)));
232
z = abs(work_matrix.coeff(p,q)) / work_matrix.coeff(p,q);
233
work_matrix.row(p) *= z;
234
if(svd.computeU()) svd.m_matrixU.col(p) *= conj(z);
235
z = abs(work_matrix.coeff(q,q)) / work_matrix.coeff(q,q);
236
work_matrix.row(q) *= z;
237
if(svd.computeU()) svd.m_matrixU.col(q) *= conj(z);
241
rot.c() = conj(work_matrix.coeff(p,p)) / n;
242
rot.s() = work_matrix.coeff(q,p) / n;
243
work_matrix.applyOnTheLeft(p,q,rot);
244
if(svd.computeU()) svd.m_matrixU.applyOnTheRight(p,q,rot.adjoint());
245
if(work_matrix.coeff(p,q) != Scalar(0))
247
Scalar z = abs(work_matrix.coeff(p,q)) / work_matrix.coeff(p,q);
248
work_matrix.col(q) *= z;
249
if(svd.computeV()) svd.m_matrixV.col(q) *= z;
251
if(work_matrix.coeff(q,q) != Scalar(0))
253
z = abs(work_matrix.coeff(q,q)) / work_matrix.coeff(q,q);
254
work_matrix.row(q) *= z;
255
if(svd.computeU()) svd.m_matrixU.col(q) *= conj(z);
261
template<typename MatrixType, typename RealScalar, typename Index>
262
void real_2x2_jacobi_svd(const MatrixType& matrix, Index p, Index q,
263
JacobiRotation<RealScalar> *j_left,
264
JacobiRotation<RealScalar> *j_right)
266
Matrix<RealScalar,2,2> m;
267
m << real(matrix.coeff(p,p)), real(matrix.coeff(p,q)),
268
real(matrix.coeff(q,p)), real(matrix.coeff(q,q));
269
JacobiRotation<RealScalar> rot1;
270
RealScalar t = m.coeff(0,0) + m.coeff(1,1);
271
RealScalar d = m.coeff(1,0) - m.coeff(0,1);
272
if(t == RealScalar(0))
274
rot1.c() = RealScalar(0);
275
rot1.s() = d > RealScalar(0) ? RealScalar(1) : RealScalar(-1);
279
RealScalar u = d / t;
280
rot1.c() = RealScalar(1) / sqrt(RealScalar(1) + abs2(u));
281
rot1.s() = rot1.c() * u;
283
m.applyOnTheLeft(0,1,rot1);
284
j_right->makeJacobi(m,0,1);
285
*j_left = rot1 * j_right->transpose();
288
} // end namespace internal
290
/** \ingroup SVD_Module
295
* \brief Two-sided Jacobi SVD decomposition of a rectangular matrix
297
* \param MatrixType the type of the matrix of which we are computing the SVD decomposition
298
* \param QRPreconditioner this optional parameter allows to specify the type of QR decomposition that will be used internally
299
* for the R-SVD step for non-square matrices. See discussion of possible values below.
301
* SVD decomposition consists in decomposing any n-by-p matrix \a A as a product
302
* \f[ A = U S V^* \f]
303
* where \a U is a n-by-n unitary, \a V is a p-by-p unitary, and \a S is a n-by-p real positive matrix which is zero outside of its main diagonal;
304
* the diagonal entries of S are known as the \em singular \em values of \a A and the columns of \a U and \a V are known as the left
305
* and right \em singular \em vectors of \a A respectively.
307
* Singular values are always sorted in decreasing order.
309
* This JacobiSVD decomposition computes only the singular values by default. If you want \a U or \a V, you need to ask for them explicitly.
311
* You can ask for only \em thin \a U or \a V to be computed, meaning the following. In case of a rectangular n-by-p matrix, letting \a m be the
312
* smaller value among \a n and \a p, there are only \a m singular vectors; the remaining columns of \a U and \a V do not correspond to actual
313
* singular vectors. Asking for \em thin \a U or \a V means asking for only their \a m first columns to be formed. So \a U is then a n-by-m matrix,
314
* and \a V is then a p-by-m matrix. Notice that thin \a U and \a V are all you need for (least squares) solving.
316
* Here's an example demonstrating basic usage:
317
* \include JacobiSVD_basic.cpp
318
* Output: \verbinclude JacobiSVD_basic.out
320
* This JacobiSVD class is a two-sided Jacobi R-SVD decomposition, ensuring optimal reliability and accuracy. The downside is that it's slower than
321
* bidiagonalizing SVD algorithms for large square matrices; however its complexity is still \f$ O(n^2p) \f$ where \a n is the smaller dimension and
322
* \a p is the greater dimension, meaning that it is still of the same order of complexity as the faster bidiagonalizing R-SVD algorithms.
323
* In particular, like any R-SVD, it takes advantage of non-squareness in that its complexity is only linear in the greater dimension.
325
* If the input matrix has inf or nan coefficients, the result of the computation is undefined, but the computation is guaranteed to
326
* terminate in finite (and reasonable) time.
328
* The possible values for QRPreconditioner are:
329
* \li ColPivHouseholderQRPreconditioner is the default. In practice it's very safe. It uses column-pivoting QR.
330
* \li FullPivHouseholderQRPreconditioner, is the safest and slowest. It uses full-pivoting QR.
331
* Contrary to other QRs, it doesn't allow computing thin unitaries.
332
* \li HouseholderQRPreconditioner is the fastest, and less safe and accurate than the pivoting variants. It uses non-pivoting QR.
333
* This is very similar in safety and accuracy to the bidiagonalization process used by bidiagonalizing SVD algorithms (since bidiagonalization
334
* is inherently non-pivoting). However the resulting SVD is still more reliable than bidiagonalizing SVDs because the Jacobi-based iterarive
335
* process is more reliable than the optimized bidiagonal SVD iterations.
336
* \li NoQRPreconditioner allows not to use a QR preconditioner at all. This is useful if you know that you will only be computing
337
* JacobiSVD decompositions of square matrices. Non-square matrices require a QR preconditioner. Using this option will result in
338
* faster compilation and smaller executable code. It won't significantly speed up computation, since JacobiSVD is always checking
339
* if QR preconditioning is needed before applying it anyway.
341
* \sa MatrixBase::jacobiSvd()
343
template<typename _MatrixType, int QRPreconditioner> class JacobiSVD
347
typedef _MatrixType MatrixType;
348
typedef typename MatrixType::Scalar Scalar;
349
typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
350
typedef typename MatrixType::Index Index;
352
RowsAtCompileTime = MatrixType::RowsAtCompileTime,
353
ColsAtCompileTime = MatrixType::ColsAtCompileTime,
354
DiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_DYNAMIC(RowsAtCompileTime,ColsAtCompileTime),
355
MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
356
MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
357
MaxDiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(MaxRowsAtCompileTime,MaxColsAtCompileTime),
358
MatrixOptions = MatrixType::Options
361
typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime,
362
MatrixOptions, MaxRowsAtCompileTime, MaxRowsAtCompileTime>
364
typedef Matrix<Scalar, ColsAtCompileTime, ColsAtCompileTime,
365
MatrixOptions, MaxColsAtCompileTime, MaxColsAtCompileTime>
367
typedef typename internal::plain_diag_type<MatrixType, RealScalar>::type SingularValuesType;
368
typedef typename internal::plain_row_type<MatrixType>::type RowType;
369
typedef typename internal::plain_col_type<MatrixType>::type ColType;
370
typedef Matrix<Scalar, DiagSizeAtCompileTime, DiagSizeAtCompileTime,
371
MatrixOptions, MaxDiagSizeAtCompileTime, MaxDiagSizeAtCompileTime>
374
/** \brief Default Constructor.
376
* The default constructor is useful in cases in which the user intends to
377
* perform decompositions via JacobiSVD::compute(const MatrixType&).
380
: m_isInitialized(false),
381
m_isAllocated(false),
382
m_computationOptions(0),
383
m_rows(-1), m_cols(-1)
387
/** \brief Default Constructor with memory preallocation
389
* Like the default constructor but with preallocation of the internal data
390
* according to the specified problem size.
393
JacobiSVD(Index rows, Index cols, unsigned int computationOptions = 0)
394
: m_isInitialized(false),
395
m_isAllocated(false),
396
m_computationOptions(0),
397
m_rows(-1), m_cols(-1)
399
allocate(rows, cols, computationOptions);
402
/** \brief Constructor performing the decomposition of given matrix.
404
* \param matrix the matrix to decompose
405
* \param computationOptions optional parameter allowing to specify if you want full or thin U or V unitaries to be computed.
406
* By default, none is computed. This is a bit-field, the possible bits are #ComputeFullU, #ComputeThinU,
407
* #ComputeFullV, #ComputeThinV.
409
* Thin unitaries are only available if your matrix type has a Dynamic number of columns (for example MatrixXf). They also are not
410
* available with the (non-default) FullPivHouseholderQR preconditioner.
412
JacobiSVD(const MatrixType& matrix, unsigned int computationOptions = 0)
413
: m_isInitialized(false),
414
m_isAllocated(false),
415
m_computationOptions(0),
416
m_rows(-1), m_cols(-1)
418
compute(matrix, computationOptions);
421
/** \brief Method performing the decomposition of given matrix using custom options.
423
* \param matrix the matrix to decompose
424
* \param computationOptions optional parameter allowing to specify if you want full or thin U or V unitaries to be computed.
425
* By default, none is computed. This is a bit-field, the possible bits are #ComputeFullU, #ComputeThinU,
426
* #ComputeFullV, #ComputeThinV.
428
* Thin unitaries are only available if your matrix type has a Dynamic number of columns (for example MatrixXf). They also are not
429
* available with the (non-default) FullPivHouseholderQR preconditioner.
431
JacobiSVD& compute(const MatrixType& matrix, unsigned int computationOptions);
433
/** \brief Method performing the decomposition of given matrix using current options.
435
* \param matrix the matrix to decompose
437
* This method uses the current \a computationOptions, as already passed to the constructor or to compute(const MatrixType&, unsigned int).
439
JacobiSVD& compute(const MatrixType& matrix)
441
return compute(matrix, m_computationOptions);
444
/** \returns the \a U matrix.
446
* For the SVD decomposition of a n-by-p matrix, letting \a m be the minimum of \a n and \a p,
447
* the U matrix is n-by-n if you asked for #ComputeFullU, and is n-by-m if you asked for #ComputeThinU.
449
* The \a m first columns of \a U are the left singular vectors of the matrix being decomposed.
451
* This method asserts that you asked for \a U to be computed.
453
const MatrixUType& matrixU() const
455
eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
456
eigen_assert(computeU() && "This JacobiSVD decomposition didn't compute U. Did you ask for it?");
460
/** \returns the \a V matrix.
462
* For the SVD decomposition of a n-by-p matrix, letting \a m be the minimum of \a n and \a p,
463
* the V matrix is p-by-p if you asked for #ComputeFullV, and is p-by-m if you asked for ComputeThinV.
465
* The \a m first columns of \a V are the right singular vectors of the matrix being decomposed.
467
* This method asserts that you asked for \a V to be computed.
469
const MatrixVType& matrixV() const
471
eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
472
eigen_assert(computeV() && "This JacobiSVD decomposition didn't compute V. Did you ask for it?");
476
/** \returns the vector of singular values.
478
* For the SVD decomposition of a n-by-p matrix, letting \a m be the minimum of \a n and \a p, the
479
* returned vector has size \a m. Singular values are always sorted in decreasing order.
481
const SingularValuesType& singularValues() const
483
eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
484
return m_singularValues;
487
/** \returns true if \a U (full or thin) is asked for in this SVD decomposition */
488
inline bool computeU() const { return m_computeFullU || m_computeThinU; }
489
/** \returns true if \a V (full or thin) is asked for in this SVD decomposition */
490
inline bool computeV() const { return m_computeFullV || m_computeThinV; }
492
/** \returns a (least squares) solution of \f$ A x = b \f$ using the current SVD decomposition of A.
494
* \param b the right-hand-side of the equation to solve.
496
* \note Solving requires both U and V to be computed. Thin U and V are enough, there is no need for full U or V.
498
* \note SVD solving is implicitly least-squares. Thus, this method serves both purposes of exact solving and least-squares solving.
499
* In other words, the returned solution is guaranteed to minimize the Euclidean norm \f$ \Vert A x - b \Vert \f$.
501
template<typename Rhs>
502
inline const internal::solve_retval<JacobiSVD, Rhs>
503
solve(const MatrixBase<Rhs>& b) const
505
eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
506
eigen_assert(computeU() && computeV() && "JacobiSVD::solve() requires both unitaries U and V to be computed (thin unitaries suffice).");
507
return internal::solve_retval<JacobiSVD, Rhs>(*this, b.derived());
510
/** \returns the number of singular values that are not exactly 0 */
511
Index nonzeroSingularValues() const
513
eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
514
return m_nonzeroSingularValues;
517
inline Index rows() const { return m_rows; }
518
inline Index cols() const { return m_cols; }
521
void allocate(Index rows, Index cols, unsigned int computationOptions);
524
MatrixUType m_matrixU;
525
MatrixVType m_matrixV;
526
SingularValuesType m_singularValues;
527
WorkMatrixType m_workMatrix;
528
bool m_isInitialized, m_isAllocated;
529
bool m_computeFullU, m_computeThinU;
530
bool m_computeFullV, m_computeThinV;
531
unsigned int m_computationOptions;
532
Index m_nonzeroSingularValues, m_rows, m_cols, m_diagSize;
534
template<typename __MatrixType, int _QRPreconditioner, bool _IsComplex>
535
friend struct internal::svd_precondition_2x2_block_to_be_real;
536
template<typename __MatrixType, int _QRPreconditioner, int _Case, bool _DoAnything>
537
friend struct internal::qr_preconditioner_impl;
540
template<typename MatrixType, int QRPreconditioner>
541
void JacobiSVD<MatrixType, QRPreconditioner>::allocate(Index rows, Index cols, unsigned int computationOptions)
543
eigen_assert(rows >= 0 && cols >= 0);
548
computationOptions == m_computationOptions)
555
m_isInitialized = false;
556
m_isAllocated = true;
557
m_computationOptions = computationOptions;
558
m_computeFullU = (computationOptions & ComputeFullU) != 0;
559
m_computeThinU = (computationOptions & ComputeThinU) != 0;
560
m_computeFullV = (computationOptions & ComputeFullV) != 0;
561
m_computeThinV = (computationOptions & ComputeThinV) != 0;
562
eigen_assert(!(m_computeFullU && m_computeThinU) && "JacobiSVD: you can't ask for both full and thin U");
563
eigen_assert(!(m_computeFullV && m_computeThinV) && "JacobiSVD: you can't ask for both full and thin V");
564
eigen_assert(EIGEN_IMPLIES(m_computeThinU || m_computeThinV, MatrixType::ColsAtCompileTime==Dynamic) &&
565
"JacobiSVD: thin U and V are only available when your matrix has a dynamic number of columns.");
566
if (QRPreconditioner == FullPivHouseholderQRPreconditioner)
568
eigen_assert(!(m_computeThinU || m_computeThinV) &&
569
"JacobiSVD: can't compute thin U or thin V with the FullPivHouseholderQR preconditioner. "
570
"Use the ColPivHouseholderQR preconditioner instead.");
572
m_diagSize = std::min(m_rows, m_cols);
573
m_singularValues.resize(m_diagSize);
574
m_matrixU.resize(m_rows, m_computeFullU ? m_rows
575
: m_computeThinU ? m_diagSize
577
m_matrixV.resize(m_cols, m_computeFullV ? m_cols
578
: m_computeThinV ? m_diagSize
580
m_workMatrix.resize(m_diagSize, m_diagSize);
583
template<typename MatrixType, int QRPreconditioner>
584
JacobiSVD<MatrixType, QRPreconditioner>&
585
JacobiSVD<MatrixType, QRPreconditioner>::compute(const MatrixType& matrix, unsigned int computationOptions)
587
allocate(matrix.rows(), matrix.cols(), computationOptions);
589
// currently we stop when we reach precision 2*epsilon as the last bit of precision can require an unreasonable number of iterations,
590
// only worsening the precision of U and V as we accumulate more rotations
591
const RealScalar precision = RealScalar(2) * NumTraits<Scalar>::epsilon();
593
/*** step 1. The R-SVD step: we use a QR decomposition to reduce to the case of a square matrix */
595
if(!internal::qr_preconditioner_impl<MatrixType, QRPreconditioner, internal::PreconditionIfMoreColsThanRows>::run(*this, matrix)
596
&& !internal::qr_preconditioner_impl<MatrixType, QRPreconditioner, internal::PreconditionIfMoreRowsThanCols>::run(*this, matrix))
598
m_workMatrix = matrix.block(0,0,m_diagSize,m_diagSize);
599
if(m_computeFullU) m_matrixU.setIdentity(m_rows,m_rows);
600
if(m_computeThinU) m_matrixU.setIdentity(m_rows,m_diagSize);
601
if(m_computeFullV) m_matrixV.setIdentity(m_cols,m_cols);
602
if(m_computeThinV) m_matrixV.setIdentity(m_cols, m_diagSize);
605
/*** step 2. The main Jacobi SVD iteration. ***/
607
bool finished = false;
612
// do a sweep: for all index pairs (p,q), perform SVD of the corresponding 2x2 sub-matrix
614
for(Index p = 1; p < m_diagSize; ++p)
616
for(Index q = 0; q < p; ++q)
618
// if this 2x2 sub-matrix is not diagonal already...
619
// notice that this comparison will evaluate to false if any NaN is involved, ensuring that NaN's don't
620
// keep us iterating forever.
622
if(max(internal::abs(m_workMatrix.coeff(p,q)),internal::abs(m_workMatrix.coeff(q,p)))
623
> max(internal::abs(m_workMatrix.coeff(p,p)),internal::abs(m_workMatrix.coeff(q,q)))*precision)
627
// perform SVD decomposition of 2x2 sub-matrix corresponding to indices p,q to make it diagonal
628
internal::svd_precondition_2x2_block_to_be_real<MatrixType, QRPreconditioner>::run(m_workMatrix, *this, p, q);
629
JacobiRotation<RealScalar> j_left, j_right;
630
internal::real_2x2_jacobi_svd(m_workMatrix, p, q, &j_left, &j_right);
632
// accumulate resulting Jacobi rotations
633
m_workMatrix.applyOnTheLeft(p,q,j_left);
634
if(computeU()) m_matrixU.applyOnTheRight(p,q,j_left.transpose());
636
m_workMatrix.applyOnTheRight(p,q,j_right);
637
if(computeV()) m_matrixV.applyOnTheRight(p,q,j_right);
643
/*** step 3. The work matrix is now diagonal, so ensure it's positive so its diagonal entries are the singular values ***/
645
for(Index i = 0; i < m_diagSize; ++i)
647
RealScalar a = internal::abs(m_workMatrix.coeff(i,i));
648
m_singularValues.coeffRef(i) = a;
649
if(computeU() && (a!=RealScalar(0))) m_matrixU.col(i) *= m_workMatrix.coeff(i,i)/a;
652
/*** step 4. Sort singular values in descending order and compute the number of nonzero singular values ***/
654
m_nonzeroSingularValues = m_diagSize;
655
for(Index i = 0; i < m_diagSize; i++)
658
RealScalar maxRemainingSingularValue = m_singularValues.tail(m_diagSize-i).maxCoeff(&pos);
659
if(maxRemainingSingularValue == RealScalar(0))
661
m_nonzeroSingularValues = i;
667
std::swap(m_singularValues.coeffRef(i), m_singularValues.coeffRef(pos));
668
if(computeU()) m_matrixU.col(pos).swap(m_matrixU.col(i));
669
if(computeV()) m_matrixV.col(pos).swap(m_matrixV.col(i));
673
m_isInitialized = true;
678
template<typename _MatrixType, int QRPreconditioner, typename Rhs>
679
struct solve_retval<JacobiSVD<_MatrixType, QRPreconditioner>, Rhs>
680
: solve_retval_base<JacobiSVD<_MatrixType, QRPreconditioner>, Rhs>
682
typedef JacobiSVD<_MatrixType, QRPreconditioner> JacobiSVDType;
683
EIGEN_MAKE_SOLVE_HELPERS(JacobiSVDType,Rhs)
685
template<typename Dest> void evalTo(Dest& dst) const
687
eigen_assert(rhs().rows() == dec().rows());
690
// So A^{-1} = V S^{-1} U^*
692
Index diagSize = std::min(dec().rows(), dec().cols());
693
typename JacobiSVDType::SingularValuesType invertedSingVals(diagSize);
695
Index nonzeroSingVals = dec().nonzeroSingularValues();
696
invertedSingVals.head(nonzeroSingVals) = dec().singularValues().head(nonzeroSingVals).array().inverse();
697
invertedSingVals.tail(diagSize - nonzeroSingVals).setZero();
699
dst = dec().matrixV().leftCols(diagSize)
700
* invertedSingVals.asDiagonal()
701
* dec().matrixU().leftCols(diagSize).adjoint()
705
} // end namespace internal
707
template<typename Derived>
708
JacobiSVD<typename MatrixBase<Derived>::PlainObject>
709
MatrixBase<Derived>::jacobiSvd(unsigned int computationOptions) const
711
return JacobiSVD<PlainObject>(*this, computationOptions);
716
#endif // EIGEN_JACOBISVD_H