Eigen  3.3.0
 
Loading...
Searching...
No Matches
JacobiSVD.h
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// Copyright (C) 2013-2014 Gael Guennebaud <gael.guennebaud@inria.fr>
6//
7// This Source Code Form is subject to the terms of the Mozilla
8// Public License v. 2.0. If a copy of the MPL was not distributed
9// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10
11#ifndef EIGEN_JACOBISVD_H
12#define EIGEN_JACOBISVD_H
13
14namespace Eigen {
15
16namespace internal {
17// forward declaration (needed by ICC)
18// the empty body is required by MSVC
19template<typename MatrixType, int QRPreconditioner,
20 bool IsComplex = NumTraits<typename MatrixType::Scalar>::IsComplex>
21struct svd_precondition_2x2_block_to_be_real {};
22
23/*** QR preconditioners (R-SVD)
24 ***
25 *** Their role is to reduce the problem of computing the SVD to the case of a square matrix.
26 *** This approach, known as R-SVD, is an optimization for rectangular-enough matrices, and is a requirement for
27 *** JacobiSVD which by itself is only able to work on square matrices.
28 ***/
29
30enum { PreconditionIfMoreColsThanRows, PreconditionIfMoreRowsThanCols };
31
32template<typename MatrixType, int QRPreconditioner, int Case>
33struct qr_preconditioner_should_do_anything
34{
35 enum { a = MatrixType::RowsAtCompileTime != Dynamic &&
36 MatrixType::ColsAtCompileTime != Dynamic &&
37 MatrixType::ColsAtCompileTime <= MatrixType::RowsAtCompileTime,
38 b = MatrixType::RowsAtCompileTime != Dynamic &&
39 MatrixType::ColsAtCompileTime != Dynamic &&
40 MatrixType::RowsAtCompileTime <= MatrixType::ColsAtCompileTime,
41 ret = !( (QRPreconditioner == NoQRPreconditioner) ||
42 (Case == PreconditionIfMoreColsThanRows && bool(a)) ||
43 (Case == PreconditionIfMoreRowsThanCols && bool(b)) )
44 };
45};
46
47template<typename MatrixType, int QRPreconditioner, int Case,
48 bool DoAnything = qr_preconditioner_should_do_anything<MatrixType, QRPreconditioner, Case>::ret
49> struct qr_preconditioner_impl {};
50
51template<typename MatrixType, int QRPreconditioner, int Case>
52class qr_preconditioner_impl<MatrixType, QRPreconditioner, Case, false>
53{
54public:
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::Scalar Scalar;
69 enum
70 {
71 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
72 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime
73 };
74 typedef Matrix<Scalar, 1, RowsAtCompileTime, RowMajor, 1, MaxRowsAtCompileTime> WorkspaceType;
75
76 void allocate(const JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner>& svd)
77 {
78 if (svd.rows() != m_qr.rows() || svd.cols() != m_qr.cols())
79 {
80 m_qr.~QRType();
81 ::new (&m_qr) QRType(svd.rows(), svd.cols());
82 }
83 if (svd.m_computeFullU) m_workspace.resize(svd.rows());
84 }
85
86 bool run(JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
87 {
88 if(matrix.rows() > matrix.cols())
89 {
90 m_qr.compute(matrix);
91 svd.m_workMatrix = m_qr.matrixQR().block(0,0,matrix.cols(),matrix.cols()).template triangularView<Upper>();
92 if(svd.m_computeFullU) m_qr.matrixQ().evalTo(svd.m_matrixU, m_workspace);
93 if(svd.computeV()) svd.m_matrixV = m_qr.colsPermutation();
94 return true;
95 }
96 return false;
97 }
98private:
99 typedef FullPivHouseholderQR<MatrixType> QRType;
100 QRType m_qr;
101 WorkspaceType m_workspace;
102};
103
104template<typename MatrixType>
105class qr_preconditioner_impl<MatrixType, FullPivHouseholderQRPreconditioner, PreconditionIfMoreColsThanRows, true>
106{
107public:
108 typedef typename MatrixType::Scalar Scalar;
109 enum
110 {
111 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
112 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
113 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
114 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
115 Options = MatrixType::Options
116 };
117 typedef Matrix<Scalar, ColsAtCompileTime, RowsAtCompileTime, Options, MaxColsAtCompileTime, MaxRowsAtCompileTime>
118 TransposeTypeWithSameStorageOrder;
119
120 void allocate(const JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner>& svd)
121 {
122 if (svd.cols() != m_qr.rows() || svd.rows() != m_qr.cols())
123 {
124 m_qr.~QRType();
125 ::new (&m_qr) QRType(svd.cols(), svd.rows());
126 }
127 m_adjoint.resize(svd.cols(), svd.rows());
128 if (svd.m_computeFullV) m_workspace.resize(svd.cols());
129 }
130
131 bool run(JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
132 {
133 if(matrix.cols() > matrix.rows())
134 {
135 m_adjoint = matrix.adjoint();
136 m_qr.compute(m_adjoint);
137 svd.m_workMatrix = m_qr.matrixQR().block(0,0,matrix.rows(),matrix.rows()).template triangularView<Upper>().adjoint();
138 if(svd.m_computeFullV) m_qr.matrixQ().evalTo(svd.m_matrixV, m_workspace);
139 if(svd.computeU()) svd.m_matrixU = m_qr.colsPermutation();
140 return true;
141 }
142 else return false;
143 }
144private:
145 typedef FullPivHouseholderQR<TransposeTypeWithSameStorageOrder> QRType;
146 QRType m_qr;
147 TransposeTypeWithSameStorageOrder m_adjoint;
148 typename internal::plain_row_type<MatrixType>::type m_workspace;
149};
150
151/*** preconditioner using ColPivHouseholderQR ***/
152
153template<typename MatrixType>
154class qr_preconditioner_impl<MatrixType, ColPivHouseholderQRPreconditioner, PreconditionIfMoreRowsThanCols, true>
155{
156public:
157 void allocate(const JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner>& svd)
158 {
159 if (svd.rows() != m_qr.rows() || svd.cols() != m_qr.cols())
160 {
161 m_qr.~QRType();
162 ::new (&m_qr) QRType(svd.rows(), svd.cols());
163 }
164 if (svd.m_computeFullU) m_workspace.resize(svd.rows());
165 else if (svd.m_computeThinU) m_workspace.resize(svd.cols());
166 }
167
168 bool run(JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
169 {
170 if(matrix.rows() > matrix.cols())
171 {
172 m_qr.compute(matrix);
173 svd.m_workMatrix = m_qr.matrixQR().block(0,0,matrix.cols(),matrix.cols()).template triangularView<Upper>();
174 if(svd.m_computeFullU) m_qr.householderQ().evalTo(svd.m_matrixU, m_workspace);
175 else if(svd.m_computeThinU)
176 {
177 svd.m_matrixU.setIdentity(matrix.rows(), matrix.cols());
178 m_qr.householderQ().applyThisOnTheLeft(svd.m_matrixU, m_workspace);
179 }
180 if(svd.computeV()) svd.m_matrixV = m_qr.colsPermutation();
181 return true;
182 }
183 return false;
184 }
185
186private:
187 typedef ColPivHouseholderQR<MatrixType> QRType;
188 QRType m_qr;
189 typename internal::plain_col_type<MatrixType>::type m_workspace;
190};
191
192template<typename MatrixType>
193class qr_preconditioner_impl<MatrixType, ColPivHouseholderQRPreconditioner, PreconditionIfMoreColsThanRows, true>
194{
195public:
196 typedef typename MatrixType::Scalar Scalar;
197 enum
198 {
199 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
200 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
201 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
202 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
203 Options = MatrixType::Options
204 };
205
206 typedef Matrix<Scalar, ColsAtCompileTime, RowsAtCompileTime, Options, MaxColsAtCompileTime, MaxRowsAtCompileTime>
207 TransposeTypeWithSameStorageOrder;
208
209 void allocate(const JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner>& svd)
210 {
211 if (svd.cols() != m_qr.rows() || svd.rows() != m_qr.cols())
212 {
213 m_qr.~QRType();
214 ::new (&m_qr) QRType(svd.cols(), svd.rows());
215 }
216 if (svd.m_computeFullV) m_workspace.resize(svd.cols());
217 else if (svd.m_computeThinV) m_workspace.resize(svd.rows());
218 m_adjoint.resize(svd.cols(), svd.rows());
219 }
220
221 bool run(JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner>& svd, const MatrixType& matrix)
222 {
223 if(matrix.cols() > matrix.rows())
224 {
225 m_adjoint = matrix.adjoint();
226 m_qr.compute(m_adjoint);
227
228 svd.m_workMatrix = m_qr.matrixQR().block(0,0,matrix.rows(),matrix.rows()).template triangularView<Upper>().adjoint();
229 if(svd.m_computeFullV) m_qr.householderQ().evalTo(svd.m_matrixV, m_workspace);
230 else if(svd.m_computeThinV)
231 {
232 svd.m_matrixV.setIdentity(matrix.cols(), matrix.rows());
233 m_qr.householderQ().applyThisOnTheLeft(svd.m_matrixV, m_workspace);
234 }
235 if(svd.computeU()) svd.m_matrixU = m_qr.colsPermutation();
236 return true;
237 }
238 else return false;
239 }
240
241private:
242 typedef ColPivHouseholderQR<TransposeTypeWithSameStorageOrder> QRType;
243 QRType m_qr;
244 TransposeTypeWithSameStorageOrder m_adjoint;
245 typename internal::plain_row_type<MatrixType>::type m_workspace;
246};
247
248/*** preconditioner using HouseholderQR ***/
249
250template<typename MatrixType>
251class qr_preconditioner_impl<MatrixType, HouseholderQRPreconditioner, PreconditionIfMoreRowsThanCols, true>
252{
253public:
254 void allocate(const JacobiSVD<MatrixType, HouseholderQRPreconditioner>& svd)
255 {
256 if (svd.rows() != m_qr.rows() || svd.cols() != m_qr.cols())
257 {
258 m_qr.~QRType();
259 ::new (&m_qr) QRType(svd.rows(), svd.cols());
260 }
261 if (svd.m_computeFullU) m_workspace.resize(svd.rows());
262 else if (svd.m_computeThinU) m_workspace.resize(svd.cols());
263 }
264
265 bool run(JacobiSVD<MatrixType, HouseholderQRPreconditioner>& svd, const MatrixType& matrix)
266 {
267 if(matrix.rows() > matrix.cols())
268 {
269 m_qr.compute(matrix);
270 svd.m_workMatrix = m_qr.matrixQR().block(0,0,matrix.cols(),matrix.cols()).template triangularView<Upper>();
271 if(svd.m_computeFullU) m_qr.householderQ().evalTo(svd.m_matrixU, m_workspace);
272 else if(svd.m_computeThinU)
273 {
274 svd.m_matrixU.setIdentity(matrix.rows(), matrix.cols());
275 m_qr.householderQ().applyThisOnTheLeft(svd.m_matrixU, m_workspace);
276 }
277 if(svd.computeV()) svd.m_matrixV.setIdentity(matrix.cols(), matrix.cols());
278 return true;
279 }
280 return false;
281 }
282private:
283 typedef HouseholderQR<MatrixType> QRType;
284 QRType m_qr;
285 typename internal::plain_col_type<MatrixType>::type m_workspace;
286};
287
288template<typename MatrixType>
289class qr_preconditioner_impl<MatrixType, HouseholderQRPreconditioner, PreconditionIfMoreColsThanRows, true>
290{
291public:
292 typedef typename MatrixType::Scalar Scalar;
293 enum
294 {
295 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
296 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
297 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
298 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
299 Options = MatrixType::Options
300 };
301
302 typedef Matrix<Scalar, ColsAtCompileTime, RowsAtCompileTime, Options, MaxColsAtCompileTime, MaxRowsAtCompileTime>
303 TransposeTypeWithSameStorageOrder;
304
305 void allocate(const JacobiSVD<MatrixType, HouseholderQRPreconditioner>& svd)
306 {
307 if (svd.cols() != m_qr.rows() || svd.rows() != m_qr.cols())
308 {
309 m_qr.~QRType();
310 ::new (&m_qr) QRType(svd.cols(), svd.rows());
311 }
312 if (svd.m_computeFullV) m_workspace.resize(svd.cols());
313 else if (svd.m_computeThinV) m_workspace.resize(svd.rows());
314 m_adjoint.resize(svd.cols(), svd.rows());
315 }
316
317 bool run(JacobiSVD<MatrixType, HouseholderQRPreconditioner>& svd, const MatrixType& matrix)
318 {
319 if(matrix.cols() > matrix.rows())
320 {
321 m_adjoint = matrix.adjoint();
322 m_qr.compute(m_adjoint);
323
324 svd.m_workMatrix = m_qr.matrixQR().block(0,0,matrix.rows(),matrix.rows()).template triangularView<Upper>().adjoint();
325 if(svd.m_computeFullV) m_qr.householderQ().evalTo(svd.m_matrixV, m_workspace);
326 else if(svd.m_computeThinV)
327 {
328 svd.m_matrixV.setIdentity(matrix.cols(), matrix.rows());
329 m_qr.householderQ().applyThisOnTheLeft(svd.m_matrixV, m_workspace);
330 }
331 if(svd.computeU()) svd.m_matrixU.setIdentity(matrix.rows(), matrix.rows());
332 return true;
333 }
334 else return false;
335 }
336
337private:
338 typedef HouseholderQR<TransposeTypeWithSameStorageOrder> QRType;
339 QRType m_qr;
340 TransposeTypeWithSameStorageOrder m_adjoint;
341 typename internal::plain_row_type<MatrixType>::type m_workspace;
342};
343
344/*** 2x2 SVD implementation
345 ***
346 *** JacobiSVD consists in performing a series of 2x2 SVD subproblems
347 ***/
348
349template<typename MatrixType, int QRPreconditioner>
350struct svd_precondition_2x2_block_to_be_real<MatrixType, QRPreconditioner, false>
351{
352 typedef JacobiSVD<MatrixType, QRPreconditioner> SVD;
353 typedef typename MatrixType::RealScalar RealScalar;
354 static bool run(typename SVD::WorkMatrixType&, SVD&, Index, Index, RealScalar&) { return true; }
355};
356
357template<typename MatrixType, int QRPreconditioner>
358struct svd_precondition_2x2_block_to_be_real<MatrixType, QRPreconditioner, true>
359{
360 typedef JacobiSVD<MatrixType, QRPreconditioner> SVD;
361 typedef typename MatrixType::Scalar Scalar;
362 typedef typename MatrixType::RealScalar RealScalar;
363 static bool run(typename SVD::WorkMatrixType& work_matrix, SVD& svd, Index p, Index q, RealScalar& maxDiagEntry)
364 {
365 using std::sqrt;
366 using std::abs;
367 Scalar z;
368 JacobiRotation<Scalar> rot;
369 RealScalar n = sqrt(numext::abs2(work_matrix.coeff(p,p)) + numext::abs2(work_matrix.coeff(q,p)));
370
371 const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
372 const RealScalar precision = NumTraits<Scalar>::epsilon();
373
374 if(n==0)
375 {
376 // make sure first column is zero
377 work_matrix.coeffRef(p,p) = work_matrix.coeffRef(q,p) = Scalar(0);
378
379 if(abs(numext::imag(work_matrix.coeff(p,q)))>considerAsZero)
380 {
381 // 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
382 z = abs(work_matrix.coeff(p,q)) / work_matrix.coeff(p,q);
383 work_matrix.row(p) *= z;
384 if(svd.computeU()) svd.m_matrixU.col(p) *= conj(z);
385 }
386 if(abs(numext::imag(work_matrix.coeff(q,q)))>considerAsZero)
387 {
388 z = abs(work_matrix.coeff(q,q)) / work_matrix.coeff(q,q);
389 work_matrix.row(q) *= z;
390 if(svd.computeU()) svd.m_matrixU.col(q) *= conj(z);
391 }
392 // otherwise the second row is already zero, so we have nothing to do.
393 }
394 else
395 {
396 rot.c() = conj(work_matrix.coeff(p,p)) / n;
397 rot.s() = work_matrix.coeff(q,p) / n;
398 work_matrix.applyOnTheLeft(p,q,rot);
399 if(svd.computeU()) svd.m_matrixU.applyOnTheRight(p,q,rot.adjoint());
400 if(abs(numext::imag(work_matrix.coeff(p,q)))>considerAsZero)
401 {
402 z = abs(work_matrix.coeff(p,q)) / work_matrix.coeff(p,q);
403 work_matrix.col(q) *= z;
404 if(svd.computeV()) svd.m_matrixV.col(q) *= z;
405 }
406 if(abs(numext::imag(work_matrix.coeff(q,q)))>considerAsZero)
407 {
408 z = abs(work_matrix.coeff(q,q)) / work_matrix.coeff(q,q);
409 work_matrix.row(q) *= z;
410 if(svd.computeU()) svd.m_matrixU.col(q) *= conj(z);
411 }
412 }
413
414 // update largest diagonal entry
415 maxDiagEntry = numext::maxi(maxDiagEntry,numext::maxi(abs(work_matrix.coeff(p,p)), abs(work_matrix.coeff(q,q))));
416 // and check whether the 2x2 block is already diagonal
417 RealScalar threshold = numext::maxi<RealScalar>(considerAsZero, precision * maxDiagEntry);
418 return abs(work_matrix.coeff(p,q))>threshold || abs(work_matrix.coeff(q,p)) > threshold;
419 }
420};
421
422template<typename _MatrixType, int QRPreconditioner>
423struct traits<JacobiSVD<_MatrixType,QRPreconditioner> >
424{
425 typedef _MatrixType MatrixType;
426};
427
428} // end namespace internal
429
483template<typename _MatrixType, int QRPreconditioner> class JacobiSVD
484 : public SVDBase<JacobiSVD<_MatrixType,QRPreconditioner> >
485{
486 typedef SVDBase<JacobiSVD> Base;
487 public:
488
489 typedef _MatrixType MatrixType;
490 typedef typename MatrixType::Scalar Scalar;
491 typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
492 enum {
493 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
494 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
495 DiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_DYNAMIC(RowsAtCompileTime,ColsAtCompileTime),
496 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
497 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
498 MaxDiagSizeAtCompileTime = EIGEN_SIZE_MIN_PREFER_FIXED(MaxRowsAtCompileTime,MaxColsAtCompileTime),
499 MatrixOptions = MatrixType::Options
500 };
501
502 typedef typename Base::MatrixUType MatrixUType;
503 typedef typename Base::MatrixVType MatrixVType;
504 typedef typename Base::SingularValuesType SingularValuesType;
505
506 typedef typename internal::plain_row_type<MatrixType>::type RowType;
507 typedef typename internal::plain_col_type<MatrixType>::type ColType;
508 typedef Matrix<Scalar, DiagSizeAtCompileTime, DiagSizeAtCompileTime,
509 MatrixOptions, MaxDiagSizeAtCompileTime, MaxDiagSizeAtCompileTime>
511
518 {}
519
520
527 JacobiSVD(Index rows, Index cols, unsigned int computationOptions = 0)
528 {
529 allocate(rows, cols, computationOptions);
530 }
531
542 explicit JacobiSVD(const MatrixType& matrix, unsigned int computationOptions = 0)
543 {
544 compute(matrix, computationOptions);
545 }
546
557 JacobiSVD& compute(const MatrixType& matrix, unsigned int computationOptions);
558
565 JacobiSVD& compute(const MatrixType& matrix)
566 {
567 return compute(matrix, m_computationOptions);
568 }
569
570 using Base::computeU;
571 using Base::computeV;
572 using Base::rows;
573 using Base::cols;
574 using Base::rank;
575
576 private:
577 void allocate(Index rows, Index cols, unsigned int computationOptions);
578
579 protected:
580 using Base::m_matrixU;
581 using Base::m_matrixV;
582 using Base::m_singularValues;
583 using Base::m_isInitialized;
584 using Base::m_isAllocated;
585 using Base::m_usePrescribedThreshold;
586 using Base::m_computeFullU;
587 using Base::m_computeThinU;
588 using Base::m_computeFullV;
589 using Base::m_computeThinV;
590 using Base::m_computationOptions;
591 using Base::m_nonzeroSingularValues;
592 using Base::m_rows;
593 using Base::m_cols;
594 using Base::m_diagSize;
595 using Base::m_prescribedThreshold;
596 WorkMatrixType m_workMatrix;
597
598 template<typename __MatrixType, int _QRPreconditioner, bool _IsComplex>
599 friend struct internal::svd_precondition_2x2_block_to_be_real;
600 template<typename __MatrixType, int _QRPreconditioner, int _Case, bool _DoAnything>
601 friend struct internal::qr_preconditioner_impl;
602
603 internal::qr_preconditioner_impl<MatrixType, QRPreconditioner, internal::PreconditionIfMoreColsThanRows> m_qr_precond_morecols;
604 internal::qr_preconditioner_impl<MatrixType, QRPreconditioner, internal::PreconditionIfMoreRowsThanCols> m_qr_precond_morerows;
605 MatrixType m_scaledMatrix;
606};
607
608template<typename MatrixType, int QRPreconditioner>
609void JacobiSVD<MatrixType, QRPreconditioner>::allocate(Index rows, Index cols, unsigned int computationOptions)
610{
611 eigen_assert(rows >= 0 && cols >= 0);
612
613 if (m_isAllocated &&
614 rows == m_rows &&
615 cols == m_cols &&
616 computationOptions == m_computationOptions)
617 {
618 return;
619 }
620
621 m_rows = rows;
622 m_cols = cols;
623 m_isInitialized = false;
624 m_isAllocated = true;
625 m_computationOptions = computationOptions;
626 m_computeFullU = (computationOptions & ComputeFullU) != 0;
627 m_computeThinU = (computationOptions & ComputeThinU) != 0;
628 m_computeFullV = (computationOptions & ComputeFullV) != 0;
629 m_computeThinV = (computationOptions & ComputeThinV) != 0;
630 eigen_assert(!(m_computeFullU && m_computeThinU) && "JacobiSVD: you can't ask for both full and thin U");
631 eigen_assert(!(m_computeFullV && m_computeThinV) && "JacobiSVD: you can't ask for both full and thin V");
632 eigen_assert(EIGEN_IMPLIES(m_computeThinU || m_computeThinV, MatrixType::ColsAtCompileTime==Dynamic) &&
633 "JacobiSVD: thin U and V are only available when your matrix has a dynamic number of columns.");
634 if (QRPreconditioner == FullPivHouseholderQRPreconditioner)
635 {
636 eigen_assert(!(m_computeThinU || m_computeThinV) &&
637 "JacobiSVD: can't compute thin U or thin V with the FullPivHouseholderQR preconditioner. "
638 "Use the ColPivHouseholderQR preconditioner instead.");
639 }
640 m_diagSize = (std::min)(m_rows, m_cols);
641 m_singularValues.resize(m_diagSize);
642 if(RowsAtCompileTime==Dynamic)
643 m_matrixU.resize(m_rows, m_computeFullU ? m_rows
644 : m_computeThinU ? m_diagSize
645 : 0);
646 if(ColsAtCompileTime==Dynamic)
647 m_matrixV.resize(m_cols, m_computeFullV ? m_cols
648 : m_computeThinV ? m_diagSize
649 : 0);
650 m_workMatrix.resize(m_diagSize, m_diagSize);
651
652 if(m_cols>m_rows) m_qr_precond_morecols.allocate(*this);
653 if(m_rows>m_cols) m_qr_precond_morerows.allocate(*this);
654 if(m_rows!=m_cols) m_scaledMatrix.resize(rows,cols);
655}
656
657template<typename MatrixType, int QRPreconditioner>
658JacobiSVD<MatrixType, QRPreconditioner>&
659JacobiSVD<MatrixType, QRPreconditioner>::compute(const MatrixType& matrix, unsigned int computationOptions)
660{
661 using std::abs;
662 allocate(matrix.rows(), matrix.cols(), computationOptions);
663
664 // currently we stop when we reach precision 2*epsilon as the last bit of precision can require an unreasonable number of iterations,
665 // only worsening the precision of U and V as we accumulate more rotations
666 const RealScalar precision = RealScalar(2) * NumTraits<Scalar>::epsilon();
667
668 // limit for denormal numbers to be considered zero in order to avoid infinite loops (see bug 286)
669 const RealScalar considerAsZero = (std::numeric_limits<RealScalar>::min)();
670
671 // Scaling factor to reduce over/under-flows
672 RealScalar scale = matrix.cwiseAbs().maxCoeff();
673 if(scale==RealScalar(0)) scale = RealScalar(1);
674
675 /*** step 1. The R-SVD step: we use a QR decomposition to reduce to the case of a square matrix */
676
677 if(m_rows!=m_cols)
678 {
679 m_scaledMatrix = matrix / scale;
680 m_qr_precond_morecols.run(*this, m_scaledMatrix);
681 m_qr_precond_morerows.run(*this, m_scaledMatrix);
682 }
683 else
684 {
685 m_workMatrix = matrix.block(0,0,m_diagSize,m_diagSize) / scale;
686 if(m_computeFullU) m_matrixU.setIdentity(m_rows,m_rows);
687 if(m_computeThinU) m_matrixU.setIdentity(m_rows,m_diagSize);
688 if(m_computeFullV) m_matrixV.setIdentity(m_cols,m_cols);
689 if(m_computeThinV) m_matrixV.setIdentity(m_cols, m_diagSize);
690 }
691
692 /*** step 2. The main Jacobi SVD iteration. ***/
693 RealScalar maxDiagEntry = m_workMatrix.cwiseAbs().diagonal().maxCoeff();
694
695 bool finished = false;
696 while(!finished)
697 {
698 finished = true;
699
700 // do a sweep: for all index pairs (p,q), perform SVD of the corresponding 2x2 sub-matrix
701
702 for(Index p = 1; p < m_diagSize; ++p)
703 {
704 for(Index q = 0; q < p; ++q)
705 {
706 // if this 2x2 sub-matrix is not diagonal already...
707 // notice that this comparison will evaluate to false if any NaN is involved, ensuring that NaN's don't
708 // keep us iterating forever. Similarly, small denormal numbers are considered zero.
709 RealScalar threshold = numext::maxi<RealScalar>(considerAsZero, precision * maxDiagEntry);
710 if(abs(m_workMatrix.coeff(p,q))>threshold || abs(m_workMatrix.coeff(q,p)) > threshold)
711 {
712 finished = false;
713 // perform SVD decomposition of 2x2 sub-matrix corresponding to indices p,q to make it diagonal
714 // the complex to real operation returns true if the updated 2x2 block is not already diagonal
715 if(internal::svd_precondition_2x2_block_to_be_real<MatrixType, QRPreconditioner>::run(m_workMatrix, *this, p, q, maxDiagEntry))
716 {
717 JacobiRotation<RealScalar> j_left, j_right;
718 internal::real_2x2_jacobi_svd(m_workMatrix, p, q, &j_left, &j_right);
719
720 // accumulate resulting Jacobi rotations
721 m_workMatrix.applyOnTheLeft(p,q,j_left);
722 if(computeU()) m_matrixU.applyOnTheRight(p,q,j_left.transpose());
723
724 m_workMatrix.applyOnTheRight(p,q,j_right);
725 if(computeV()) m_matrixV.applyOnTheRight(p,q,j_right);
726
727 // keep track of the largest diagonal coefficient
728 maxDiagEntry = numext::maxi<RealScalar>(maxDiagEntry,numext::maxi(abs(m_workMatrix.coeff(p,p)), abs(m_workMatrix.coeff(q,q))));
729 }
730 }
731 }
732 }
733 }
734
735 /*** step 3. The work matrix is now diagonal, so ensure it's positive so its diagonal entries are the singular values ***/
736
737 for(Index i = 0; i < m_diagSize; ++i)
738 {
739 // For a complex matrix, some diagonal coefficients might note have been
740 // treated by svd_precondition_2x2_block_to_be_real, and the imaginary part
741 // of some diagonal entry might not be null.
742 if(NumTraits<Scalar>::IsComplex && abs(numext::imag(m_workMatrix.coeff(i,i)))>considerAsZero)
743 {
744 RealScalar a = abs(m_workMatrix.coeff(i,i));
745 m_singularValues.coeffRef(i) = abs(a);
746 if(computeU()) m_matrixU.col(i) *= m_workMatrix.coeff(i,i)/a;
747 }
748 else
749 {
750 // m_workMatrix.coeff(i,i) is already real, no difficulty:
751 RealScalar a = numext::real(m_workMatrix.coeff(i,i));
752 m_singularValues.coeffRef(i) = abs(a);
753 if(computeU() && (a<RealScalar(0))) m_matrixU.col(i) = -m_matrixU.col(i);
754 }
755 }
756
757 m_singularValues *= scale;
758
759 /*** step 4. Sort singular values in descending order and compute the number of nonzero singular values ***/
760
761 m_nonzeroSingularValues = m_diagSize;
762 for(Index i = 0; i < m_diagSize; i++)
763 {
764 Index pos;
765 RealScalar maxRemainingSingularValue = m_singularValues.tail(m_diagSize-i).maxCoeff(&pos);
766 if(maxRemainingSingularValue == RealScalar(0))
767 {
768 m_nonzeroSingularValues = i;
769 break;
770 }
771 if(pos)
772 {
773 pos += i;
774 std::swap(m_singularValues.coeffRef(i), m_singularValues.coeffRef(pos));
775 if(computeU()) m_matrixU.col(pos).swap(m_matrixU.col(i));
776 if(computeV()) m_matrixV.col(pos).swap(m_matrixV.col(i));
777 }
778 }
779
780 m_isInitialized = true;
781 return *this;
782}
783
791template<typename Derived>
793MatrixBase<Derived>::jacobiSvd(unsigned int computationOptions) const
794{
795 return JacobiSVD<PlainObject>(*this, computationOptions);
796}
797
798} // end namespace Eigen
799
800#endif // EIGEN_JACOBISVD_H
internal::traits< Derived >::Scalar Scalar
Definition: DenseBase.h:66
Rotation given by a cosine-sine pair.
Definition: Jacobi.h:35
JacobiRotation transpose() const
Definition: Jacobi.h:59
Two-sided Jacobi SVD decomposition of a rectangular matrix.
Definition: JacobiSVD.h:485
JacobiSVD & compute(const MatrixType &matrix)
Method performing the decomposition of given matrix using current options.
Definition: JacobiSVD.h:565
JacobiSVD()
Default Constructor.
Definition: JacobiSVD.h:517
JacobiSVD(Index rows, Index cols, unsigned int computationOptions=0)
Default Constructor with memory preallocation.
Definition: JacobiSVD.h:527
JacobiSVD & compute(const MatrixType &matrix, unsigned int computationOptions)
Method performing the decomposition of given matrix using custom options.
Definition: JacobiSVD.h:659
JacobiSVD(const MatrixType &matrix, unsigned int computationOptions=0)
Constructor performing the decomposition of given matrix.
Definition: JacobiSVD.h:542
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:50
Base class of SVD algorithms.
Definition: SVDBase.h:49
Index rank() const
Definition: SVDBase.h:130
bool computeV() const
Definition: SVDBase.h:190
bool computeU() const
Definition: SVDBase.h:188
@ NoQRPreconditioner
Definition: Constants.h:415
@ HouseholderQRPreconditioner
Definition: Constants.h:417
@ ColPivHouseholderQRPreconditioner
Definition: Constants.h:419
@ FullPivHouseholderQRPreconditioner
Definition: Constants.h:421
@ ComputeFullV
Definition: Constants.h:387
@ ComputeThinV
Definition: Constants.h:389
@ ComputeFullU
Definition: Constants.h:383
@ ComputeThinU
Definition: Constants.h:385
Namespace containing all symbols from the Eigen library.
Definition: Core:287
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_conjugate_op< typename Derived::Scalar >, const Derived > conj(const Eigen::ArrayBase< Derived > &x)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
const int Dynamic
Definition: Constants.h:21
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:151