source: pacpussensors/trunk/Vislab/lib3dv/eigen/Eigen/src/SVD/JacobiSVD.h@ 136

Last change on this file since 136 was 136, checked in by ldecherf, 7 years ago

Doc

File size: 39.9 KB
Line 
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2009-2010 Benoit Jacob <jacob.benoit.1@gmail.com>
5//
6// This Source Code Form is subject to the terms of the Mozilla
7// Public License v. 2.0. If a copy of the MPL was not distributed
8// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9
10#ifndef EIGEN_JACOBISVD_H
11#define EIGEN_JACOBISVD_H
12
13namespace Eigen {
14
15namespace internal {
16// forward declaration (needed by ICC)
17// the empty body is required by MSVC
18template<typename MatrixType, int QRPreconditioner,
19 bool IsComplex = NumTraits<typename MatrixType::Scalar>::IsComplex>
20struct svd_precondition_2x2_block_to_be_real {};
21
22/*** QR preconditioners (R-SVD)
23 ***
24 *** Their role is to reduce the problem of computing the SVD to the case of a square matrix.
25 *** This approach, known as R-SVD, is an optimization for rectangular-enough matrices, and is a requirement for
26 *** JacobiSVD which by itself is only able to work on square matrices.
27 ***/
28
29enum { PreconditionIfMoreColsThanRows, PreconditionIfMoreRowsThanCols };
30
31template<typename MatrixType, int QRPreconditioner, int Case>
32struct qr_preconditioner_should_do_anything
33{
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)) )
43 };
44};
45
46template<typename MatrixType, int QRPreconditioner, int Case,
47 bool DoAnything = qr_preconditioner_should_do_anything<MatrixType, QRPreconditioner, Case>::ret
48> struct qr_preconditioner_impl {};
49
50template<typename MatrixType, int QRPreconditioner, int Case>
51class qr_preconditioner_impl<MatrixType, QRPreconditioner, Case, false>
52{
53public:
54 typedef typename MatrixType::Index Index;
55 void allocate(const JacobiSVD<MatrixType, QRPreconditioner>&) {}
56 bool run(JacobiSVD<MatrixType, QRPreconditioner>&, const MatrixType&)
57 {
58 return false;
59 }
60};
61
62/*** preconditioner using FullPivHouseholderQR ***/
63
64template<typename MatrixType>
65class qr_preconditioner_impl<MatrixType, FullPivHouseholderQRPreconditioner, PreconditionIfMoreRowsThanCols, true>
66{
67public:
68 typedef typename MatrixType::Index Index;
69 typedef typename MatrixType::Scalar Scalar;
70 enum
71 {
72 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
73 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime
74 };
75 typedef Matrix<Scalar, 1, RowsAtCompileTime, RowMajor, 1, MaxRowsAtCompileTime> WorkspaceType;
76
77 void allocate(const JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner>& svd)
78 {
79 if (svd.rows() != m_qr.rows() || svd.cols() != m_qr.cols())
80 {
81 m_qr.~QRType();
82 ::new (&m_qr) QRType(svd.rows(), svd.cols());
83 }
84 if (svd.m_computeFullU) m_workspace.resize(svd.rows());
85 }
86
87 bool run(JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
88 {
89 if(matrix.rows() > matrix.cols())
90 {
91 m_qr.compute(matrix);
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();
95 return true;
96 }
97 return false;
98 }
99private:
100 typedef FullPivHouseholderQR<MatrixType> QRType;
101 QRType m_qr;
102 WorkspaceType m_workspace;
103};
104
105template<typename MatrixType>
106class qr_preconditioner_impl<MatrixType, FullPivHouseholderQRPreconditioner, PreconditionIfMoreColsThanRows, true>
107{
108public:
109 typedef typename MatrixType::Index Index;
110 typedef typename MatrixType::Scalar Scalar;
111 enum
112 {
113 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
114 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
115 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
116 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
117 Options = MatrixType::Options
118 };
119 typedef Matrix<Scalar, ColsAtCompileTime, RowsAtCompileTime, Options, MaxColsAtCompileTime, MaxRowsAtCompileTime>
120 TransposeTypeWithSameStorageOrder;
121
122 void allocate(const JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner>& svd)
123 {
124 if (svd.cols() != m_qr.rows() || svd.rows() != m_qr.cols())
125 {
126 m_qr.~QRType();
127 ::new (&m_qr) QRType(svd.cols(), svd.rows());
128 }
129 m_adjoint.resize(svd.cols(), svd.rows());
130 if (svd.m_computeFullV) m_workspace.resize(svd.cols());
131 }
132
133 bool run(JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
134 {
135 if(matrix.cols() > matrix.rows())
136 {
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();
142 return true;
143 }
144 else return false;
145 }
146private:
147 typedef FullPivHouseholderQR<TransposeTypeWithSameStorageOrder> QRType;
148 QRType m_qr;
149 TransposeTypeWithSameStorageOrder m_adjoint;
150 typename internal::plain_row_type<MatrixType>::type m_workspace;
151};
152
153/*** preconditioner using ColPivHouseholderQR ***/
154
155template<typename MatrixType>
156class qr_preconditioner_impl<MatrixType, ColPivHouseholderQRPreconditioner, PreconditionIfMoreRowsThanCols, true>
157{
158public:
159 typedef typename MatrixType::Index Index;
160
161 void allocate(const JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner>& svd)
162 {
163 if (svd.rows() != m_qr.rows() || svd.cols() != m_qr.cols())
164 {
165 m_qr.~QRType();
166 ::new (&m_qr) QRType(svd.rows(), svd.cols());
167 }
168 if (svd.m_computeFullU) m_workspace.resize(svd.rows());
169 else if (svd.m_computeThinU) m_workspace.resize(svd.cols());
170 }
171
172 bool run(JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
173 {
174 if(matrix.rows() > matrix.cols())
175 {
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)
180 {
181 svd.m_matrixU.setIdentity(matrix.rows(), matrix.cols());
182 m_qr.householderQ().applyThisOnTheLeft(svd.m_matrixU, m_workspace);
183 }
184 if(svd.computeV()) svd.m_matrixV = m_qr.colsPermutation();
185 return true;
186 }
187 return false;
188 }
189
190private:
191 typedef ColPivHouseholderQR<MatrixType> QRType;
192 QRType m_qr;
193 typename internal::plain_col_type<MatrixType>::type m_workspace;
194};
195
196template<typename MatrixType>
197class qr_preconditioner_impl<MatrixType, ColPivHouseholderQRPreconditioner, PreconditionIfMoreColsThanRows, true>
198{
199public:
200 typedef typename MatrixType::Index Index;
201 typedef typename MatrixType::Scalar Scalar;
202 enum
203 {
204 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
205 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
206 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
207 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
208 Options = MatrixType::Options
209 };
210
211 typedef Matrix<Scalar, ColsAtCompileTime, RowsAtCompileTime, Options, MaxColsAtCompileTime, MaxRowsAtCompileTime>
212 TransposeTypeWithSameStorageOrder;
213
214 void allocate(const JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner>& svd)
215 {
216 if (svd.cols() != m_qr.rows() || svd.rows() != m_qr.cols())
217 {
218 m_qr.~QRType();
219 ::new (&m_qr) QRType(svd.cols(), svd.rows());
220 }
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());
224 }
225
226 bool run(JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
227 {
228 if(matrix.cols() > matrix.rows())
229 {
230 m_adjoint = matrix.adjoint();
231 m_qr.compute(m_adjoint);
232
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)
236 {
237 svd.m_matrixV.setIdentity(matrix.cols(), matrix.rows());
238 m_qr.householderQ().applyThisOnTheLeft(svd.m_matrixV, m_workspace);
239 }
240 if(svd.computeU()) svd.m_matrixU = m_qr.colsPermutation();
241 return true;
242 }
243 else return false;
244 }
245
246private:
247 typedef ColPivHouseholderQR<TransposeTypeWithSameStorageOrder> QRType;
248 QRType m_qr;
249 TransposeTypeWithSameStorageOrder m_adjoint;
250 typename internal::plain_row_type<MatrixType>::type m_workspace;
251};
252
253/*** preconditioner using HouseholderQR ***/
254
255template<typename MatrixType>
256class qr_preconditioner_impl<MatrixType, HouseholderQRPreconditioner, PreconditionIfMoreRowsThanCols, true>
257{
258public:
259 typedef typename MatrixType::Index Index;
260
261 void allocate(const JacobiSVD<MatrixType, HouseholderQRPreconditioner>& svd)
262 {
263 if (svd.rows() != m_qr.rows() || svd.cols() != m_qr.cols())
264 {
265 m_qr.~QRType();
266 ::new (&m_qr) QRType(svd.rows(), svd.cols());
267 }
268 if (svd.m_computeFullU) m_workspace.resize(svd.rows());
269 else if (svd.m_computeThinU) m_workspace.resize(svd.cols());
270 }
271
272 bool run(JacobiSVD<MatrixType, HouseholderQRPreconditioner>& svd, const MatrixType& matrix)
273 {
274 if(matrix.rows() > matrix.cols())
275 {
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)
280 {
281 svd.m_matrixU.setIdentity(matrix.rows(), matrix.cols());
282 m_qr.householderQ().applyThisOnTheLeft(svd.m_matrixU, m_workspace);
283 }
284 if(svd.computeV()) svd.m_matrixV.setIdentity(matrix.cols(), matrix.cols());
285 return true;
286 }
287 return false;
288 }
289private:
290 typedef HouseholderQR<MatrixType> QRType;
291 QRType m_qr;
292 typename internal::plain_col_type<MatrixType>::type m_workspace;
293};
294
295template<typename MatrixType>
296class qr_preconditioner_impl<MatrixType, HouseholderQRPreconditioner, PreconditionIfMoreColsThanRows, true>
297{
298public:
299 typedef typename MatrixType::Index Index;
300 typedef typename MatrixType::Scalar Scalar;
301 enum
302 {
303 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
304 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
305 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
306 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
307 Options = MatrixType::Options
308 };
309
310 typedef Matrix<Scalar, ColsAtCompileTime, RowsAtCompileTime, Options, MaxColsAtCompileTime, MaxRowsAtCompileTime>
311 TransposeTypeWithSameStorageOrder;
312
313 void allocate(const JacobiSVD<MatrixType, HouseholderQRPreconditioner>& svd)
314 {
315 if (svd.cols() != m_qr.rows() || svd.rows() != m_qr.cols())
316 {
317 m_qr.~QRType();
318 ::new (&m_qr) QRType(svd.cols(), svd.rows());
319 }
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());
323 }
324
325 bool run(JacobiSVD<MatrixType, HouseholderQRPreconditioner>& svd, const MatrixType& matrix)
326 {
327 if(matrix.cols() > matrix.rows())
328 {
329 m_adjoint = matrix.adjoint();
330 m_qr.compute(m_adjoint);
331
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)
335 {
336 svd.m_matrixV.setIdentity(matrix.cols(), matrix.rows());
337 m_qr.householderQ().applyThisOnTheLeft(svd.m_matrixV, m_workspace);
338 }
339 if(svd.computeU()) svd.m_matrixU.setIdentity(matrix.rows(), matrix.rows());
340 return true;
341 }
342 else return false;
343 }
344
345private:
346 typedef HouseholderQR<TransposeTypeWithSameStorageOrder> QRType;
347 QRType m_qr;
348 TransposeTypeWithSameStorageOrder m_adjoint;
349 typename internal::plain_row_type<MatrixType>::type m_workspace;
350};
351
352/*** 2x2 SVD implementation
353 ***
354 *** JacobiSVD consists in performing a series of 2x2 SVD subproblems
355 ***/
356
357template<typename MatrixType, int QRPreconditioner>
358struct svd_precondition_2x2_block_to_be_real<MatrixType, QRPreconditioner, false>
359{
360 typedef JacobiSVD<MatrixType, QRPreconditioner> SVD;
361 typedef typename SVD::Index Index;
362 typedef typename MatrixType::RealScalar RealScalar;
363 static bool run(typename SVD::WorkMatrixType&, SVD&, Index, Index, RealScalar&) { return true; }
364};
365
366template<typename MatrixType, int QRPreconditioner>
367struct svd_precondition_2x2_block_to_be_real<MatrixType, QRPreconditioner, true>
368{
369 typedef JacobiSVD<MatrixType, QRPreconditioner> SVD;
370 typedef typename SVD::Index Index;
371 typedef typename MatrixType::Scalar Scalar;
372 typedef typename MatrixType::RealScalar RealScalar;
373 static bool run(typename SVD::WorkMatrixType& work_matrix, SVD& svd, Index p, Index q, RealScalar& maxDiagEntry)
374 {
375 using std::sqrt;
376 using std::abs;
377 using std::max;
378 Scalar z;
379 JacobiRotation<Scalar> rot;
380 RealScalar n = sqrt(numext::abs2(work_matrix.coeff(p,p)) + numext::abs2(work_matrix.coeff(q,p)));
381
382 const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
383 const RealScalar precision = NumTraits<Scalar>::epsilon();
384
385 if(n==0)
386 {
387 // make sure first column is zero
388 work_matrix.coeffRef(p,p) = work_matrix.coeffRef(q,p) = Scalar(0);
389
390 if(abs(numext::imag(work_matrix.coeff(p,q)))>considerAsZero)
391 {
392 // work_matrix.coeff(p,q) can be zero if work_matrix.coeff(q,p) is not zero but small enough to underflow when computing n
393 z = abs(work_matrix.coeff(p,q)) / work_matrix.coeff(p,q);
394 work_matrix.row(p) *= z;
395 if(svd.computeU()) svd.m_matrixU.col(p) *= conj(z);
396 }
397 if(abs(numext::imag(work_matrix.coeff(q,q)))>considerAsZero)
398 {
399 z = abs(work_matrix.coeff(q,q)) / work_matrix.coeff(q,q);
400 work_matrix.row(q) *= z;
401 if(svd.computeU()) svd.m_matrixU.col(q) *= conj(z);
402 }
403 // otherwise the second row is already zero, so we have nothing to do.
404 }
405 else
406 {
407 rot.c() = conj(work_matrix.coeff(p,p)) / n;
408 rot.s() = work_matrix.coeff(q,p) / n;
409 work_matrix.applyOnTheLeft(p,q,rot);
410 if(svd.computeU()) svd.m_matrixU.applyOnTheRight(p,q,rot.adjoint());
411 if(abs(numext::imag(work_matrix.coeff(p,q)))>considerAsZero)
412 {
413 z = abs(work_matrix.coeff(p,q)) / work_matrix.coeff(p,q);
414 work_matrix.col(q) *= z;
415 if(svd.computeV()) svd.m_matrixV.col(q) *= z;
416 }
417 if(abs(numext::imag(work_matrix.coeff(q,q)))>considerAsZero)
418 {
419 z = abs(work_matrix.coeff(q,q)) / work_matrix.coeff(q,q);
420 work_matrix.row(q) *= z;
421 if(svd.computeU()) svd.m_matrixU.col(q) *= conj(z);
422 }
423 }
424
425 // update largest diagonal entry
426 maxDiagEntry = max EIGEN_EMPTY (maxDiagEntry,max EIGEN_EMPTY (abs(work_matrix.coeff(p,p)), abs(work_matrix.coeff(q,q))));
427 // and check whether the 2x2 block is already diagonal
428 RealScalar threshold = max EIGEN_EMPTY (considerAsZero, precision * maxDiagEntry);
429 return abs(work_matrix.coeff(p,q))>threshold || abs(work_matrix.coeff(q,p)) > threshold;
430 }
431};
432
433template<typename MatrixType, typename RealScalar, typename Index>
434void real_2x2_jacobi_svd(const MatrixType& matrix, Index p, Index q,
435 JacobiRotation<RealScalar> *j_left,
436 JacobiRotation<RealScalar> *j_right)
437{
438 using std::sqrt;
439 using std::abs;
440 Matrix<RealScalar,2,2> m;
441 m << numext::real(matrix.coeff(p,p)), numext::real(matrix.coeff(p,q)),
442 numext::real(matrix.coeff(q,p)), numext::real(matrix.coeff(q,q));
443 JacobiRotation<RealScalar> rot1;
444 RealScalar t = m.coeff(0,0) + m.coeff(1,1);
445 RealScalar d = m.coeff(1,0) - m.coeff(0,1);
446 if(d == RealScalar(0))
447 {
448 rot1.s() = RealScalar(0);
449 rot1.c() = RealScalar(1);
450 }
451 else
452 {
453 // If d!=0, then t/d cannot overflow because the magnitude of the
454 // entries forming d are not too small compared to the ones forming t.
455 RealScalar u = t / d;
456 RealScalar tmp = sqrt(RealScalar(1) + numext::abs2(u));
457 rot1.s() = RealScalar(1) / tmp;
458 rot1.c() = u / tmp;
459 }
460 m.applyOnTheLeft(0,1,rot1);
461 j_right->makeJacobi(m,0,1);
462 *j_left = rot1 * j_right->transpose();
463}
464
465} // end namespace internal
466
467/** \ingroup SVD_Module
468 *
469 *
470 * \class JacobiSVD
471 *
472 * \brief Two-sided Jacobi SVD decomposition of a rectangular matrix
473 *
474 * \param MatrixType the type of the matrix of which we are computing the SVD decomposition
475 * \param QRPreconditioner this optional parameter allows to specify the type of QR decomposition that will be used internally
476 * for the R-SVD step for non-square matrices. See discussion of possible values below.
477 *
478 * SVD decomposition consists in decomposing any n-by-p matrix \a A as a product
479 * \f[ A = U S V^* \f]
480 * 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;
481 * 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
482 * and right \em singular \em vectors of \a A respectively.
483 *
484 * Singular values are always sorted in decreasing order.
485 *
486 * 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.
487 *
488 * 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
489 * 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
490 * 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,
491 * 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.
492 *
493 * Here's an example demonstrating basic usage:
494 * \include JacobiSVD_basic.cpp
495 * Output: \verbinclude JacobiSVD_basic.out
496 *
497 * This JacobiSVD class is a two-sided Jacobi R-SVD decomposition, ensuring optimal reliability and accuracy. The downside is that it's slower than
498 * 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
499 * \a p is the greater dimension, meaning that it is still of the same order of complexity as the faster bidiagonalizing R-SVD algorithms.
500 * In particular, like any R-SVD, it takes advantage of non-squareness in that its complexity is only linear in the greater dimension.
501 *
502 * If the input matrix has inf or nan coefficients, the result of the computation is undefined, but the computation is guaranteed to
503 * terminate in finite (and reasonable) time.
504 *
505 * The possible values for QRPreconditioner are:
506 * \li ColPivHouseholderQRPreconditioner is the default. In practice it's very safe. It uses column-pivoting QR.
507 * \li FullPivHouseholderQRPreconditioner, is the safest and slowest. It uses full-pivoting QR.
508 * Contrary to other QRs, it doesn't allow computing thin unitaries.
509 * \li HouseholderQRPreconditioner is the fastest, and less safe and accurate than the pivoting variants. It uses non-pivoting QR.
510 * This is very similar in safety and accuracy to the bidiagonalization process used by bidiagonalizing SVD algorithms (since bidiagonalization
511 * is inherently non-pivoting). However the resulting SVD is still more reliable than bidiagonalizing SVDs because the Jacobi-based iterarive
512 * process is more reliable than the optimized bidiagonal SVD iterations.
513 * \li NoQRPreconditioner allows not to use a QR preconditioner at all. This is useful if you know that you will only be computing
514 * JacobiSVD decompositions of square matrices. Non-square matrices require a QR preconditioner. Using this option will result in
515 * faster compilation and smaller executable code. It won't significantly speed up computation, since JacobiSVD is always checking
516 * if QR preconditioning is needed before applying it anyway.
517 *
518 * \sa MatrixBase::jacobiSvd()
519 */
520template<typename _MatrixType, int QRPreconditioner> class JacobiSVD
521{
522 public:
523
524 typedef _MatrixType MatrixType;
525 typedef typename MatrixType::Scalar Scalar;
526 typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
527 typedef typename MatrixType::Index Index;
528 enum {
529 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
530 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
531 DiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_DYNAMIC(RowsAtCompileTime,ColsAtCompileTime),
532 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
533 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
534 MaxDiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(MaxRowsAtCompileTime,MaxColsAtCompileTime),
535 MatrixOptions = MatrixType::Options
536 };
537
538 typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime,
539 MatrixOptions, MaxRowsAtCompileTime, MaxRowsAtCompileTime>
540 MatrixUType;
541 typedef Matrix<Scalar, ColsAtCompileTime, ColsAtCompileTime,
542 MatrixOptions, MaxColsAtCompileTime, MaxColsAtCompileTime>
543 MatrixVType;
544 typedef typename internal::plain_diag_type<MatrixType, RealScalar>::type SingularValuesType;
545 typedef typename internal::plain_row_type<MatrixType>::type RowType;
546 typedef typename internal::plain_col_type<MatrixType>::type ColType;
547 typedef Matrix<Scalar, DiagSizeAtCompileTime, DiagSizeAtCompileTime,
548 MatrixOptions, MaxDiagSizeAtCompileTime, MaxDiagSizeAtCompileTime>
549 WorkMatrixType;
550
551 /** \brief Default Constructor.
552 *
553 * The default constructor is useful in cases in which the user intends to
554 * perform decompositions via JacobiSVD::compute(const MatrixType&).
555 */
556 JacobiSVD()
557 : m_isInitialized(false),
558 m_isAllocated(false),
559 m_usePrescribedThreshold(false),
560 m_computationOptions(0),
561 m_rows(-1), m_cols(-1), m_diagSize(0)
562 {}
563
564
565 /** \brief Default Constructor with memory preallocation
566 *
567 * Like the default constructor but with preallocation of the internal data
568 * according to the specified problem size.
569 * \sa JacobiSVD()
570 */
571 JacobiSVD(Index rows, Index cols, unsigned int computationOptions = 0)
572 : m_isInitialized(false),
573 m_isAllocated(false),
574 m_usePrescribedThreshold(false),
575 m_computationOptions(0),
576 m_rows(-1), m_cols(-1)
577 {
578 allocate(rows, cols, computationOptions);
579 }
580
581 /** \brief Constructor performing the decomposition of given matrix.
582 *
583 * \param matrix the matrix to decompose
584 * \param computationOptions optional parameter allowing to specify if you want full or thin U or V unitaries to be computed.
585 * By default, none is computed. This is a bit-field, the possible bits are #ComputeFullU, #ComputeThinU,
586 * #ComputeFullV, #ComputeThinV.
587 *
588 * Thin unitaries are only available if your matrix type has a Dynamic number of columns (for example MatrixXf). They also are not
589 * available with the (non-default) FullPivHouseholderQR preconditioner.
590 */
591 JacobiSVD(const MatrixType& matrix, unsigned int computationOptions = 0)
592 : m_isInitialized(false),
593 m_isAllocated(false),
594 m_usePrescribedThreshold(false),
595 m_computationOptions(0),
596 m_rows(-1), m_cols(-1)
597 {
598 compute(matrix, computationOptions);
599 }
600
601 /** \brief Method performing the decomposition of given matrix using custom options.
602 *
603 * \param matrix the matrix to decompose
604 * \param computationOptions optional parameter allowing to specify if you want full or thin U or V unitaries to be computed.
605 * By default, none is computed. This is a bit-field, the possible bits are #ComputeFullU, #ComputeThinU,
606 * #ComputeFullV, #ComputeThinV.
607 *
608 * Thin unitaries are only available if your matrix type has a Dynamic number of columns (for example MatrixXf). They also are not
609 * available with the (non-default) FullPivHouseholderQR preconditioner.
610 */
611 JacobiSVD& compute(const MatrixType& matrix, unsigned int computationOptions);
612
613 /** \brief Method performing the decomposition of given matrix using current options.
614 *
615 * \param matrix the matrix to decompose
616 *
617 * This method uses the current \a computationOptions, as already passed to the constructor or to compute(const MatrixType&, unsigned int).
618 */
619 JacobiSVD& compute(const MatrixType& matrix)
620 {
621 return compute(matrix, m_computationOptions);
622 }
623
624 /** \returns the \a U matrix.
625 *
626 * For the SVD decomposition of a n-by-p matrix, letting \a m be the minimum of \a n and \a p,
627 * the U matrix is n-by-n if you asked for #ComputeFullU, and is n-by-m if you asked for #ComputeThinU.
628 *
629 * The \a m first columns of \a U are the left singular vectors of the matrix being decomposed.
630 *
631 * This method asserts that you asked for \a U to be computed.
632 */
633 const MatrixUType& matrixU() const
634 {
635 eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
636 eigen_assert(computeU() && "This JacobiSVD decomposition didn't compute U. Did you ask for it?");
637 return m_matrixU;
638 }
639
640 /** \returns the \a V matrix.
641 *
642 * For the SVD decomposition of a n-by-p matrix, letting \a m be the minimum of \a n and \a p,
643 * the V matrix is p-by-p if you asked for #ComputeFullV, and is p-by-m if you asked for ComputeThinV.
644 *
645 * The \a m first columns of \a V are the right singular vectors of the matrix being decomposed.
646 *
647 * This method asserts that you asked for \a V to be computed.
648 */
649 const MatrixVType& matrixV() const
650 {
651 eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
652 eigen_assert(computeV() && "This JacobiSVD decomposition didn't compute V. Did you ask for it?");
653 return m_matrixV;
654 }
655
656 /** \returns the vector of singular values.
657 *
658 * For the SVD decomposition of a n-by-p matrix, letting \a m be the minimum of \a n and \a p, the
659 * returned vector has size \a m. Singular values are always sorted in decreasing order.
660 */
661 const SingularValuesType& singularValues() const
662 {
663 eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
664 return m_singularValues;
665 }
666
667 /** \returns true if \a U (full or thin) is asked for in this SVD decomposition */
668 inline bool computeU() const { return m_computeFullU || m_computeThinU; }
669 /** \returns true if \a V (full or thin) is asked for in this SVD decomposition */
670 inline bool computeV() const { return m_computeFullV || m_computeThinV; }
671
672 /** \returns a (least squares) solution of \f$ A x = b \f$ using the current SVD decomposition of A.
673 *
674 * \param b the right-hand-side of the equation to solve.
675 *
676 * \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.
677 *
678 * \note SVD solving is implicitly least-squares. Thus, this method serves both purposes of exact solving and least-squares solving.
679 * In other words, the returned solution is guaranteed to minimize the Euclidean norm \f$ \Vert A x - b \Vert \f$.
680 */
681 template<typename Rhs>
682 inline const internal::solve_retval<JacobiSVD, Rhs>
683 solve(const MatrixBase<Rhs>& b) const
684 {
685 eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
686 eigen_assert(computeU() && computeV() && "JacobiSVD::solve() requires both unitaries U and V to be computed (thin unitaries suffice).");
687 return internal::solve_retval<JacobiSVD, Rhs>(*this, b.derived());
688 }
689
690 /** \returns the number of singular values that are not exactly 0 */
691 Index nonzeroSingularValues() const
692 {
693 eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
694 return m_nonzeroSingularValues;
695 }
696
697 /** \returns the rank of the matrix of which \c *this is the SVD.
698 *
699 * \note This method has to determine which singular values should be considered nonzero.
700 * For that, it uses the threshold value that you can control by calling
701 * setThreshold(const RealScalar&).
702 */
703 inline Index rank() const
704 {
705 using std::abs;
706 eigen_assert(m_isInitialized && "JacobiSVD is not initialized.");
707 if(m_singularValues.size()==0) return 0;
708 RealScalar premultiplied_threshold = m_singularValues.coeff(0) * threshold();
709 Index i = m_nonzeroSingularValues-1;
710 while(i>=0 && m_singularValues.coeff(i) < premultiplied_threshold) --i;
711 return i+1;
712 }
713
714 /** Allows to prescribe a threshold to be used by certain methods, such as rank() and solve(),
715 * which need to determine when singular values are to be considered nonzero.
716 * This is not used for the SVD decomposition itself.
717 *
718 * When it needs to get the threshold value, Eigen calls threshold().
719 * The default is \c NumTraits<Scalar>::epsilon()
720 *
721 * \param threshold The new value to use as the threshold.
722 *
723 * A singular value will be considered nonzero if its value is strictly greater than
724 * \f$ \vert singular value \vert \leqslant threshold \times \vert max singular value \vert \f$.
725 *
726 * If you want to come back to the default behavior, call setThreshold(Default_t)
727 */
728 JacobiSVD& setThreshold(const RealScalar& threshold)
729 {
730 m_usePrescribedThreshold = true;
731 m_prescribedThreshold = threshold;
732 return *this;
733 }
734
735 /** Allows to come back to the default behavior, letting Eigen use its default formula for
736 * determining the threshold.
737 *
738 * You should pass the special object Eigen::Default as parameter here.
739 * \code svd.setThreshold(Eigen::Default); \endcode
740 *
741 * See the documentation of setThreshold(const RealScalar&).
742 */
743 JacobiSVD& setThreshold(Default_t)
744 {
745 m_usePrescribedThreshold = false;
746 return *this;
747 }
748
749 /** Returns the threshold that will be used by certain methods such as rank().
750 *
751 * See the documentation of setThreshold(const RealScalar&).
752 */
753 RealScalar threshold() const
754 {
755 eigen_assert(m_isInitialized || m_usePrescribedThreshold);
756 return m_usePrescribedThreshold ? m_prescribedThreshold
757 : (std::max<Index>)(1,m_diagSize)*NumTraits<Scalar>::epsilon();
758 }
759
760 inline Index rows() const { return m_rows; }
761 inline Index cols() const { return m_cols; }
762
763 private:
764 void allocate(Index rows, Index cols, unsigned int computationOptions);
765
766 static void check_template_parameters()
767 {
768 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
769 }
770
771 protected:
772 MatrixUType m_matrixU;
773 MatrixVType m_matrixV;
774 SingularValuesType m_singularValues;
775 WorkMatrixType m_workMatrix;
776 bool m_isInitialized, m_isAllocated, m_usePrescribedThreshold;
777 bool m_computeFullU, m_computeThinU;
778 bool m_computeFullV, m_computeThinV;
779 unsigned int m_computationOptions;
780 Index m_nonzeroSingularValues, m_rows, m_cols, m_diagSize;
781 RealScalar m_prescribedThreshold;
782
783 template<typename __MatrixType, int _QRPreconditioner, bool _IsComplex>
784 friend struct internal::svd_precondition_2x2_block_to_be_real;
785 template<typename __MatrixType, int _QRPreconditioner, int _Case, bool _DoAnything>
786 friend struct internal::qr_preconditioner_impl;
787
788 internal::qr_preconditioner_impl<MatrixType, QRPreconditioner, internal::PreconditionIfMoreColsThanRows> m_qr_precond_morecols;
789 internal::qr_preconditioner_impl<MatrixType, QRPreconditioner, internal::PreconditionIfMoreRowsThanCols> m_qr_precond_morerows;
790 MatrixType m_scaledMatrix;
791};
792
793template<typename MatrixType, int QRPreconditioner>
794void JacobiSVD<MatrixType, QRPreconditioner>::allocate(Index rows, Index cols, unsigned int computationOptions)
795{
796 eigen_assert(rows >= 0 && cols >= 0);
797
798 if (m_isAllocated &&
799 rows == m_rows &&
800 cols == m_cols &&
801 computationOptions == m_computationOptions)
802 {
803 return;
804 }
805
806 m_rows = rows;
807 m_cols = cols;
808 m_isInitialized = false;
809 m_isAllocated = true;
810 m_computationOptions = computationOptions;
811 m_computeFullU = (computationOptions & ComputeFullU) != 0;
812 m_computeThinU = (computationOptions & ComputeThinU) != 0;
813 m_computeFullV = (computationOptions & ComputeFullV) != 0;
814 m_computeThinV = (computationOptions & ComputeThinV) != 0;
815 eigen_assert(!(m_computeFullU && m_computeThinU) && "JacobiSVD: you can't ask for both full and thin U");
816 eigen_assert(!(m_computeFullV && m_computeThinV) && "JacobiSVD: you can't ask for both full and thin V");
817 eigen_assert(EIGEN_IMPLIES(m_computeThinU || m_computeThinV, MatrixType::ColsAtCompileTime==Dynamic) &&
818 "JacobiSVD: thin U and V are only available when your matrix has a dynamic number of columns.");
819 if (QRPreconditioner == FullPivHouseholderQRPreconditioner)
820 {
821 eigen_assert(!(m_computeThinU || m_computeThinV) &&
822 "JacobiSVD: can't compute thin U or thin V with the FullPivHouseholderQR preconditioner. "
823 "Use the ColPivHouseholderQR preconditioner instead.");
824 }
825 m_diagSize = (std::min)(m_rows, m_cols);
826 m_singularValues.resize(m_diagSize);
827 if(RowsAtCompileTime==Dynamic)
828 m_matrixU.resize(m_rows, m_computeFullU ? m_rows
829 : m_computeThinU ? m_diagSize
830 : 0);
831 if(ColsAtCompileTime==Dynamic)
832 m_matrixV.resize(m_cols, m_computeFullV ? m_cols
833 : m_computeThinV ? m_diagSize
834 : 0);
835 m_workMatrix.resize(m_diagSize, m_diagSize);
836
837 if(m_cols>m_rows) m_qr_precond_morecols.allocate(*this);
838 if(m_rows>m_cols) m_qr_precond_morerows.allocate(*this);
839 if(m_rows!=m_cols) m_scaledMatrix.resize(rows,cols);
840}
841
842template<typename MatrixType, int QRPreconditioner>
843JacobiSVD<MatrixType, QRPreconditioner>&
844JacobiSVD<MatrixType, QRPreconditioner>::compute(const MatrixType& matrix, unsigned int computationOptions)
845{
846 check_template_parameters();
847
848 using std::abs;
849 using std::max;
850 allocate(matrix.rows(), matrix.cols(), computationOptions);
851
852 // currently we stop when we reach precision 2*epsilon as the last bit of precision can require an unreasonable number of iterations,
853 // only worsening the precision of U and V as we accumulate more rotations
854 const RealScalar precision = RealScalar(2) * NumTraits<Scalar>::epsilon();
855
856 // limit for very small denormal numbers to be considered zero in order to avoid infinite loops (see bug 286)
857 const RealScalar considerAsZero = RealScalar(2) * std::numeric_limits<RealScalar>::denorm_min();
858
859 // Scaling factor to reduce over/under-flows
860 RealScalar scale = matrix.cwiseAbs().maxCoeff();
861 if(scale==RealScalar(0)) scale = RealScalar(1);
862
863 /*** step 1. The R-SVD step: we use a QR decomposition to reduce to the case of a square matrix */
864
865 if(m_rows!=m_cols)
866 {
867 m_scaledMatrix = matrix / scale;
868 m_qr_precond_morecols.run(*this, m_scaledMatrix);
869 m_qr_precond_morerows.run(*this, m_scaledMatrix);
870 }
871 else
872 {
873 m_workMatrix = matrix.block(0,0,m_diagSize,m_diagSize) / scale;
874 if(m_computeFullU) m_matrixU.setIdentity(m_rows,m_rows);
875 if(m_computeThinU) m_matrixU.setIdentity(m_rows,m_diagSize);
876 if(m_computeFullV) m_matrixV.setIdentity(m_cols,m_cols);
877 if(m_computeThinV) m_matrixV.setIdentity(m_cols, m_diagSize);
878 }
879
880 /*** step 2. The main Jacobi SVD iteration. ***/
881 RealScalar maxDiagEntry = m_workMatrix.cwiseAbs().diagonal().maxCoeff();
882
883 bool finished = false;
884 while(!finished)
885 {
886 finished = true;
887
888 // do a sweep: for all index pairs (p,q), perform SVD of the corresponding 2x2 sub-matrix
889
890 for(Index p = 1; p < m_diagSize; ++p)
891 {
892 for(Index q = 0; q < p; ++q)
893 {
894 // if this 2x2 sub-matrix is not diagonal already...
895 // notice that this comparison will evaluate to false if any NaN is involved, ensuring that NaN's don't
896 // keep us iterating forever. Similarly, small denormal numbers are considered zero.
897 RealScalar threshold = max EIGEN_EMPTY (considerAsZero, precision * maxDiagEntry);
898 if(abs(m_workMatrix.coeff(p,q))>threshold || abs(m_workMatrix.coeff(q,p)) > threshold)
899 {
900 finished = false;
901 // perform SVD decomposition of 2x2 sub-matrix corresponding to indices p,q to make it diagonal
902 // the complex to real operation returns true is the updated 2x2 block is not already diagonal
903 if(internal::svd_precondition_2x2_block_to_be_real<MatrixType, QRPreconditioner>::run(m_workMatrix, *this, p, q, maxDiagEntry))
904 {
905 JacobiRotation<RealScalar> j_left, j_right;
906 internal::real_2x2_jacobi_svd(m_workMatrix, p, q, &j_left, &j_right);
907
908 // accumulate resulting Jacobi rotations
909 m_workMatrix.applyOnTheLeft(p,q,j_left);
910 if(computeU()) m_matrixU.applyOnTheRight(p,q,j_left.transpose());
911
912 m_workMatrix.applyOnTheRight(p,q,j_right);
913 if(computeV()) m_matrixV.applyOnTheRight(p,q,j_right);
914
915 // keep track of the largest diagonal coefficient
916 maxDiagEntry = max EIGEN_EMPTY (maxDiagEntry,max EIGEN_EMPTY (abs(m_workMatrix.coeff(p,p)), abs(m_workMatrix.coeff(q,q))));
917 }
918 }
919 }
920 }
921 }
922
923 /*** step 3. The work matrix is now diagonal, so ensure it's positive so its diagonal entries are the singular values ***/
924
925 for(Index i = 0; i < m_diagSize; ++i)
926 {
927 RealScalar a = abs(m_workMatrix.coeff(i,i));
928 m_singularValues.coeffRef(i) = a;
929 if(computeU() && (a!=RealScalar(0))) m_matrixU.col(i) *= m_workMatrix.coeff(i,i)/a;
930 }
931
932 /*** step 4. Sort singular values in descending order and compute the number of nonzero singular values ***/
933
934 m_nonzeroSingularValues = m_diagSize;
935 for(Index i = 0; i < m_diagSize; i++)
936 {
937 Index pos;
938 RealScalar maxRemainingSingularValue = m_singularValues.tail(m_diagSize-i).maxCoeff(&pos);
939 if(maxRemainingSingularValue == RealScalar(0))
940 {
941 m_nonzeroSingularValues = i;
942 break;
943 }
944 if(pos)
945 {
946 pos += i;
947 std::swap(m_singularValues.coeffRef(i), m_singularValues.coeffRef(pos));
948 if(computeU()) m_matrixU.col(pos).swap(m_matrixU.col(i));
949 if(computeV()) m_matrixV.col(pos).swap(m_matrixV.col(i));
950 }
951 }
952
953 m_singularValues *= scale;
954
955 m_isInitialized = true;
956 return *this;
957}
958
959namespace internal {
960template<typename _MatrixType, int QRPreconditioner, typename Rhs>
961struct solve_retval<JacobiSVD<_MatrixType, QRPreconditioner>, Rhs>
962 : solve_retval_base<JacobiSVD<_MatrixType, QRPreconditioner>, Rhs>
963{
964 typedef JacobiSVD<_MatrixType, QRPreconditioner> JacobiSVDType;
965 EIGEN_MAKE_SOLVE_HELPERS(JacobiSVDType,Rhs)
966
967 template<typename Dest> void evalTo(Dest& dst) const
968 {
969 eigen_assert(rhs().rows() == dec().rows());
970
971 // A = U S V^*
972 // So A^{-1} = V S^{-1} U^*
973
974 Matrix<Scalar, Dynamic, Rhs::ColsAtCompileTime, 0, _MatrixType::MaxRowsAtCompileTime, Rhs::MaxColsAtCompileTime> tmp;
975 Index rank = dec().rank();
976
977 tmp.noalias() = dec().matrixU().leftCols(rank).adjoint() * rhs();
978 tmp = dec().singularValues().head(rank).asDiagonal().inverse() * tmp;
979 dst = dec().matrixV().leftCols(rank) * tmp;
980 }
981};
982} // end namespace internal
983
984/** \svd_module
985 *
986 * \return the singular value decomposition of \c *this computed by two-sided
987 * Jacobi transformations.
988 *
989 * \sa class JacobiSVD
990 */
991template<typename Derived>
992JacobiSVD<typename MatrixBase<Derived>::PlainObject>
993MatrixBase<Derived>::jacobiSvd(unsigned int computationOptions) const
994{
995 return JacobiSVD<PlainObject>(*this, computationOptions);
996}
997
998} // end namespace Eigen
999
1000#endif // EIGEN_JACOBISVD_H
Note: See TracBrowser for help on using the repository browser.