1
// This file is part of Eigen, a lightweight C++ template library
2
// for linear algebra. Eigen itself is part of the KDE project.
4
// Copyright (C) 2008 Gael Guennebaud <g.gael@free.fr>
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/>.
28
/** \ingroup SVD_Module
33
* \brief Standard SVD decomposition of a matrix and associated features
35
* \param MatrixType the type of the matrix of which we are computing the SVD decomposition
37
* This class performs a standard SVD decomposition of a real matrix A of size \c M x \c N
41
* \sa MatrixBase::SVD()
43
template<typename MatrixType> class SVD
46
typedef typename MatrixType::Scalar Scalar;
47
typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
50
PacketSize = internal::packet_traits<Scalar>::size,
51
AlignmentMask = int(PacketSize)-1,
52
MinSize = EIGEN_SIZE_MIN_PREFER_DYNAMIC(MatrixType::RowsAtCompileTime, MatrixType::ColsAtCompileTime)
55
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> ColVector;
56
typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, 1> RowVector;
58
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MinSize> MatrixUType;
59
typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, MatrixType::ColsAtCompileTime> MatrixVType;
60
typedef Matrix<Scalar, MinSize, 1> SingularValuesType;
64
SVD() {} // a user who relied on compiler-generated default compiler reported problems with MSVC in 2.0.7
66
SVD(const MatrixType& matrix)
67
: m_matU(matrix.rows(), std::min(matrix.rows(), matrix.cols())),
68
m_matV(matrix.cols(),matrix.cols()),
69
m_sigma(std::min(matrix.rows(),matrix.cols()))
74
template<typename OtherDerived, typename ResultType>
75
bool solve(const MatrixBase<OtherDerived> &b, ResultType* result) const;
77
const MatrixUType& matrixU() const { return m_matU; }
78
const SingularValuesType& singularValues() const { return m_sigma; }
79
const MatrixVType& matrixV() const { return m_matV; }
81
void compute(const MatrixType& matrix);
84
template<typename UnitaryType, typename PositiveType>
85
void computeUnitaryPositive(UnitaryType *unitary, PositiveType *positive) const;
86
template<typename PositiveType, typename UnitaryType>
87
void computePositiveUnitary(PositiveType *positive, UnitaryType *unitary) const;
88
template<typename RotationType, typename ScalingType>
89
void computeRotationScaling(RotationType *unitary, ScalingType *positive) const;
90
template<typename ScalingType, typename RotationType>
91
void computeScalingRotation(ScalingType *positive, RotationType *unitary) const;
99
SingularValuesType m_sigma;
102
/** Computes / recomputes the SVD decomposition A = U S V^* of \a matrix
104
* \note this code has been adapted from JAMA (public domain)
106
template<typename MatrixType>
107
void SVD<MatrixType>::compute(const MatrixType& matrix)
109
const int m = matrix.rows();
110
const int n = matrix.cols();
111
const int nu = std::min(m,n);
112
ei_assert(m>=n && "In Eigen 2.0, SVD only works for MxN matrices with M>=N. Sorry!");
113
ei_assert(m>1 && "In Eigen 2.0, SVD doesn't work on 1x1 matrices");
115
m_matU.resize(m, nu);
117
m_sigma.resize(std::min(m,n));
122
MatrixType matA(matrix);
123
const bool wantu = true;
124
const bool wantv = true;
127
// Reduce A to bidiagonal form, storing the diagonal elements
128
// in s and the super-diagonal elements in e.
129
int nct = std::min(m-1,n);
130
int nrt = std::max(0,std::min(n-2,m));
131
for (k = 0; k < std::max(nct,nrt); ++k)
135
// Compute the transformation for the k-th column and
136
// place the k-th diagonal in m_sigma[k].
137
m_sigma[k] = matA.col(k).end(m-k).norm();
138
if (m_sigma[k] != 0.0) // FIXME
141
m_sigma[k] = -m_sigma[k];
142
matA.col(k).end(m-k) /= m_sigma[k];
145
m_sigma[k] = -m_sigma[k];
148
for (j = k+1; j < n; ++j)
150
if ((k < nct) && (m_sigma[k] != 0.0))
152
// Apply the transformation.
153
Scalar t = matA.col(k).end(m-k).eigen2_dot(matA.col(j).end(m-k)); // FIXME dot product or cwise prod + .sum() ??
155
matA.col(j).end(m-k) += t * matA.col(k).end(m-k);
158
// Place the k-th row of A into e for the
159
// subsequent calculation of the row transformation.
163
// Place the transformation in U for subsequent back multiplication.
164
if (wantu & (k < nct))
165
m_matU.col(k).end(m-k) = matA.col(k).end(m-k);
169
// Compute the k-th row transformation and place the
170
// k-th super-diagonal in e[k].
171
e[k] = e.end(n-k-1).norm();
176
e.end(n-k-1) /= e[k];
180
if ((k+1 < m) & (e[k] != 0.0))
182
// Apply the transformation.
183
work.end(m-k-1) = matA.corner(BottomRight,m-k-1,n-k-1) * e.end(n-k-1);
184
for (j = k+1; j < n; ++j)
185
matA.col(j).end(m-k-1) += (-e[j]/e[k+1]) * work.end(m-k-1);
188
// Place the transformation in V for subsequent back multiplication.
190
m_matV.col(k).end(n-k-1) = e.end(n-k-1);
195
// Set up the final bidiagonal matrix or order p.
196
int p = std::min(n,m+1);
198
m_sigma[nct] = matA(nct,nct);
202
e[nrt] = matA(nrt,p-1);
205
// If required, generate U.
208
for (j = nct; j < nu; ++j)
210
m_matU.col(j).setZero();
213
for (k = nct-1; k >= 0; k--)
215
if (m_sigma[k] != 0.0)
217
for (j = k+1; j < nu; ++j)
219
Scalar t = m_matU.col(k).end(m-k).eigen2_dot(m_matU.col(j).end(m-k)); // FIXME is it really a dot product we want ?
221
m_matU.col(j).end(m-k) += t * m_matU.col(k).end(m-k);
223
m_matU.col(k).end(m-k) = - m_matU.col(k).end(m-k);
224
m_matU(k,k) = Scalar(1) + m_matU(k,k);
226
m_matU.col(k).start(k-1).setZero();
230
m_matU.col(k).setZero();
236
// If required, generate V.
239
for (k = n-1; k >= 0; k--)
241
if ((k < nrt) & (e[k] != 0.0))
243
for (j = k+1; j < nu; ++j)
245
Scalar t = m_matV.col(k).end(n-k-1).eigen2_dot(m_matV.col(j).end(n-k-1)); // FIXME is it really a dot product we want ?
246
t = -t/m_matV(k+1,k);
247
m_matV.col(j).end(n-k-1) += t * m_matV.col(k).end(n-k-1);
250
m_matV.col(k).setZero();
255
// Main iteration loop for the singular values.
258
Scalar eps = ei_pow(Scalar(2),ei_is_same_type<Scalar,float>::ret ? Scalar(-23) : Scalar(-52));
264
// Here is where a test for too many iterations would go.
266
// This section of the program inspects for
267
// negligible elements in the s and e arrays. On
268
// completion the variables kase and k are set as follows.
270
// kase = 1 if s(p) and e[k-1] are negligible and k<p
271
// kase = 2 if s(k) is negligible and k<p
272
// kase = 3 if e[k-1] is negligible, k<p, and
273
// s(k), ..., s(p) are not negligible (qr step).
274
// kase = 4 if e(p-1) is negligible (convergence).
276
for (k = p-2; k >= -1; --k)
280
if (ei_abs(e[k]) <= eps*(ei_abs(m_sigma[k]) + ei_abs(m_sigma[k+1])))
293
for (ks = p-1; ks >= k; --ks)
297
Scalar t = (ks != p ? ei_abs(e[ks]) : Scalar(0)) + (ks != k+1 ? ei_abs(e[ks-1]) : Scalar(0));
298
if (ei_abs(m_sigma[ks]) <= eps*t)
320
// Perform the task indicated by kase.
324
// Deflate negligible s(p).
329
for (j = p-2; j >= k; --j)
331
Scalar t(internal::hypot(m_sigma[j],f));
332
Scalar cs(m_sigma[j]/t);
342
for (i = 0; i < n; ++i)
344
t = cs*m_matV(i,j) + sn*m_matV(i,p-1);
345
m_matV(i,p-1) = -sn*m_matV(i,j) + cs*m_matV(i,p-1);
353
// Split at negligible s(k).
358
for (j = k; j < p; ++j)
360
Scalar t(internal::hypot(m_sigma[j],f));
361
Scalar cs( m_sigma[j]/t);
368
for (i = 0; i < m; ++i)
370
t = cs*m_matU(i,j) + sn*m_matU(i,k-1);
371
m_matU(i,k-1) = -sn*m_matU(i,j) + cs*m_matU(i,k-1);
379
// Perform one qr step.
382
// Calculate the shift.
383
Scalar scale = std::max(std::max(std::max(std::max(
384
ei_abs(m_sigma[p-1]),ei_abs(m_sigma[p-2])),ei_abs(e[p-2])),
385
ei_abs(m_sigma[k])),ei_abs(e[k]));
386
Scalar sp = m_sigma[p-1]/scale;
387
Scalar spm1 = m_sigma[p-2]/scale;
388
Scalar epm1 = e[p-2]/scale;
389
Scalar sk = m_sigma[k]/scale;
390
Scalar ek = e[k]/scale;
391
Scalar b = ((spm1 + sp)*(spm1 - sp) + epm1*epm1)/Scalar(2);
392
Scalar c = (sp*epm1)*(sp*epm1);
394
if ((b != 0.0) || (c != 0.0))
396
shift = ei_sqrt(b*b + c);
399
shift = c/(b + shift);
401
Scalar f = (sk + sp)*(sk - sp) + shift;
406
for (j = k; j < p-1; ++j)
408
Scalar t = internal::hypot(f,g);
413
f = cs*m_sigma[j] + sn*e[j];
414
e[j] = cs*e[j] - sn*m_sigma[j];
416
m_sigma[j+1] = cs*m_sigma[j+1];
419
for (i = 0; i < n; ++i)
421
t = cs*m_matV(i,j) + sn*m_matV(i,j+1);
422
m_matV(i,j+1) = -sn*m_matV(i,j) + cs*m_matV(i,j+1);
426
t = internal::hypot(f,g);
430
f = cs*e[j] + sn*m_sigma[j+1];
431
m_sigma[j+1] = -sn*e[j] + cs*m_sigma[j+1];
434
if (wantu && (j < m-1))
436
for (i = 0; i < m; ++i)
438
t = cs*m_matU(i,j) + sn*m_matU(i,j+1);
439
m_matU(i,j+1) = -sn*m_matU(i,j) + cs*m_matU(i,j+1);
452
// Make the singular values positive.
453
if (m_sigma[k] <= 0.0)
455
m_sigma[k] = m_sigma[k] < Scalar(0) ? -m_sigma[k] : Scalar(0);
457
m_matV.col(k).start(pp+1) = -m_matV.col(k).start(pp+1);
460
// Order the singular values.
463
if (m_sigma[k] >= m_sigma[k+1])
465
Scalar t = m_sigma[k];
466
m_sigma[k] = m_sigma[k+1];
468
if (wantv && (k < n-1))
469
m_matV.col(k).swap(m_matV.col(k+1));
470
if (wantu && (k < m-1))
471
m_matU.col(k).swap(m_matU.col(k+1));
482
template<typename MatrixType>
483
SVD<MatrixType>& SVD<MatrixType>::sort()
485
int mu = m_matU.rows();
486
int mv = m_matV.rows();
487
int n = m_matU.cols();
489
for (int i=0; i<n; ++i)
492
Scalar p = m_sigma.coeff(i);
494
for (int j=i+1; j<n; ++j)
496
if (m_sigma.coeff(j) > p)
499
p = m_sigma.coeff(j);
504
m_sigma.coeffRef(k) = m_sigma.coeff(i); // i.e.
505
m_sigma.coeffRef(i) = p; // swaps the i-th and the k-th elements
508
for(int s=0; j!=0; ++s, --j)
509
std::swap(m_matU.coeffRef(s,i), m_matU.coeffRef(s,k));
512
for (int s=0; j!=0; ++s, --j)
513
std::swap(m_matV.coeffRef(s,i), m_matV.coeffRef(s,k));
519
/** \returns the solution of \f$ A x = b \f$ using the current SVD decomposition of A.
520
* The parts of the solution corresponding to zero singular values are ignored.
522
* \sa MatrixBase::svd(), LU::solve(), LLT::solve()
524
template<typename MatrixType>
525
template<typename OtherDerived, typename ResultType>
526
bool SVD<MatrixType>::solve(const MatrixBase<OtherDerived> &b, ResultType* result) const
528
const int rows = m_matU.rows();
529
ei_assert(b.rows() == rows);
531
Scalar maxVal = m_sigma.cwise().abs().maxCoeff();
532
for (int j=0; j<b.cols(); ++j)
534
Matrix<Scalar,MatrixUType::RowsAtCompileTime,1> aux = m_matU.transpose() * b.col(j);
536
for (int i = 0; i <m_matU.cols(); ++i)
538
Scalar si = m_sigma.coeff(i);
539
if (ei_isMuchSmallerThan(ei_abs(si),maxVal))
542
aux.coeffRef(i) /= si;
545
result->col(j) = m_matV * aux;
550
/** Computes the polar decomposition of the matrix, as a product unitary x positive.
552
* If either pointer is zero, the corresponding computation is skipped.
554
* Only for square matrices.
556
* \sa computePositiveUnitary(), computeRotationScaling()
558
template<typename MatrixType>
559
template<typename UnitaryType, typename PositiveType>
560
void SVD<MatrixType>::computeUnitaryPositive(UnitaryType *unitary,
561
PositiveType *positive) const
563
ei_assert(m_matU.cols() == m_matV.cols() && "Polar decomposition is only for square matrices");
564
if(unitary) *unitary = m_matU * m_matV.adjoint();
565
if(positive) *positive = m_matV * m_sigma.asDiagonal() * m_matV.adjoint();
568
/** Computes the polar decomposition of the matrix, as a product positive x unitary.
570
* If either pointer is zero, the corresponding computation is skipped.
572
* Only for square matrices.
574
* \sa computeUnitaryPositive(), computeRotationScaling()
576
template<typename MatrixType>
577
template<typename UnitaryType, typename PositiveType>
578
void SVD<MatrixType>::computePositiveUnitary(UnitaryType *positive,
579
PositiveType *unitary) const
581
ei_assert(m_matU.rows() == m_matV.rows() && "Polar decomposition is only for square matrices");
582
if(unitary) *unitary = m_matU * m_matV.adjoint();
583
if(positive) *positive = m_matU * m_sigma.asDiagonal() * m_matU.adjoint();
586
/** decomposes the matrix as a product rotation x scaling, the scaling being
587
* not necessarily positive.
589
* If either pointer is zero, the corresponding computation is skipped.
591
* This method requires the Geometry module.
593
* \sa computeScalingRotation(), computeUnitaryPositive()
595
template<typename MatrixType>
596
template<typename RotationType, typename ScalingType>
597
void SVD<MatrixType>::computeRotationScaling(RotationType *rotation, ScalingType *scaling) const
599
ei_assert(m_matU.rows() == m_matV.rows() && "Polar decomposition is only for square matrices");
600
Scalar x = (m_matU * m_matV.adjoint()).determinant(); // so x has absolute value 1
601
Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> sv(m_sigma);
603
if(scaling) scaling->lazyAssign(m_matV * sv.asDiagonal() * m_matV.adjoint());
606
MatrixType m(m_matU);
608
rotation->lazyAssign(m * m_matV.adjoint());
612
/** decomposes the matrix as a product scaling x rotation, the scaling being
613
* not necessarily positive.
615
* If either pointer is zero, the corresponding computation is skipped.
617
* This method requires the Geometry module.
619
* \sa computeRotationScaling(), computeUnitaryPositive()
621
template<typename MatrixType>
622
template<typename ScalingType, typename RotationType>
623
void SVD<MatrixType>::computeScalingRotation(ScalingType *scaling, RotationType *rotation) const
625
ei_assert(m_matU.rows() == m_matV.rows() && "Polar decomposition is only for square matrices");
626
Scalar x = (m_matU * m_matV.adjoint()).determinant(); // so x has absolute value 1
627
Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> sv(m_sigma);
629
if(scaling) scaling->lazyAssign(m_matU * sv.asDiagonal() * m_matU.adjoint());
632
MatrixType m(m_matU);
634
rotation->lazyAssign(m * m_matV.adjoint());
640
* \returns the SVD decomposition of \c *this
642
template<typename Derived>
643
inline SVD<typename MatrixBase<Derived>::PlainObject>
644
MatrixBase<Derived>::svd() const
646
return SVD<PlainObject>(derived());
649
#endif // EIGEN2_SVD_H