1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
5 // Copyright (C) 2010,2012 Jitse Niesen <jitse@maths.leeds.ac.uk>
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_EIGENSOLVER_H
12 #define EIGEN_EIGENSOLVER_H
13
14 #include "./RealSchur.h"
15
16 namespace Eigen {
17
18 /** \eigenvalues_module \ingroup Eigenvalues_Module
19 *
20 *
21 * \class EigenSolver
22 *
23 * \brief Computes eigenvalues and eigenvectors of general matrices
24 *
25 * \tparam _MatrixType the type of the matrix of which we are computing the
26 * eigendecomposition; this is expected to be an instantiation of the Matrix
27 * class template. Currently, only real matrices are supported.
28 *
29 * The eigenvalues and eigenvectors of a matrix \f$ A \f$ are scalars
30 * \f$ \lambda \f$ and vectors \f$ v \f$ such that \f$ Av = \lambda v \f$. If
31 * \f$ D \f$ is a diagonal matrix with the eigenvalues on the diagonal, and
32 * \f$ V \f$ is a matrix with the eigenvectors as its columns, then \f$ A V =
33 * V D \f$. The matrix \f$ V \f$ is almost always invertible, in which case we
34 * have \f$ A = V D V^{-1} \f$. This is called the eigendecomposition.
35 *
36 * The eigenvalues and eigenvectors of a matrix may be complex, even when the
37 * matrix is real. However, we can choose real matrices \f$ V \f$ and \f$ D
38 * \f$ satisfying \f$ A V = V D \f$, just like the eigendecomposition, if the
39 * matrix \f$ D \f$ is not required to be diagonal, but if it is allowed to
40 * have blocks of the form
41 * \f[ \begin{bmatrix} u & v \\ -v & u \end{bmatrix} \f]
42 * (where \f$ u \f$ and \f$ v \f$ are real numbers) on the diagonal. These
43 * blocks correspond to complex eigenvalue pairs \f$ u \pm iv \f$. We call
44 * this variant of the eigendecomposition the pseudo-eigendecomposition.
45 *
46 * Call the function compute() to compute the eigenvalues and eigenvectors of
47 * a given matrix. Alternatively, you can use the
48 * EigenSolver(const MatrixType&, bool) constructor which computes the
49 * eigenvalues and eigenvectors at construction time. Once the eigenvalue and
50 * eigenvectors are computed, they can be retrieved with the eigenvalues() and
51 * eigenvectors() functions. The pseudoEigenvalueMatrix() and
52 * pseudoEigenvectors() methods allow the construction of the
53 * pseudo-eigendecomposition.
54 *
55 * The documentation for EigenSolver(const MatrixType&, bool) contains an
56 * example of the typical use of this class.
57 *
58 * \note The implementation is adapted from
59 * <a href="http://math.nist.gov/javanumerics/jama/">JAMA</a> (public domain).
60 * Their code is based on EISPACK.
61 *
62 * \sa MatrixBase::eigenvalues(), class ComplexEigenSolver, class SelfAdjointEigenSolver
63 */
64 template<typename _MatrixType> class EigenSolver
65 {
66 public:
67
68 /** \brief Synonym for the template parameter \p _MatrixType. */
69 typedef _MatrixType MatrixType;
70
71 enum {
72 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
73 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
74 Options = MatrixType::Options,
75 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
76 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
77 };
78
79 /** \brief Scalar type for matrices of type #MatrixType. */
80 typedef typename MatrixType::Scalar Scalar;
81 typedef typename NumTraits<Scalar>::Real RealScalar;
82 typedef Eigen::Index Index; ///< \deprecated since Eigen 3.3
83
84 /** \brief Complex scalar type for #MatrixType.
85 *
86 * This is \c std::complex<Scalar> if #Scalar is real (e.g.,
87 * \c float or \c double) and just \c Scalar if #Scalar is
88 * complex.
89 */
90 typedef std::complex<RealScalar> ComplexScalar;
91
92 /** \brief Type for vector of eigenvalues as returned by eigenvalues().
93 *
94 * This is a column vector with entries of type #ComplexScalar.
95 * The length of the vector is the size of #MatrixType.
96 */
97 typedef Matrix<ComplexScalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> EigenvalueType;
98
99 /** \brief Type for matrix of eigenvectors as returned by eigenvectors().
100 *
101 * This is a square matrix with entries of type #ComplexScalar.
102 * The size is the same as the size of #MatrixType.
103 */
104 typedef Matrix<ComplexScalar, RowsAtCompileTime, ColsAtCompileTime, Options, MaxRowsAtCompileTime, MaxColsAtCompileTime> EigenvectorsType;
105
106 /** \brief Default constructor.
107 *
108 * The default constructor is useful in cases in which the user intends to
109 * perform decompositions via EigenSolver::compute(const MatrixType&, bool).
110 *
111 * \sa compute() for an example.
112 */
EigenSolver()113 EigenSolver() : m_eivec(), m_eivalues(), m_isInitialized(false), m_realSchur(), m_matT(), m_tmp() {}
114
115 /** \brief Default constructor with memory preallocation
116 *
117 * Like the default constructor but with preallocation of the internal data
118 * according to the specified problem \a size.
119 * \sa EigenSolver()
120 */
EigenSolver(Index size)121 explicit EigenSolver(Index size)
122 : m_eivec(size, size),
123 m_eivalues(size),
124 m_isInitialized(false),
125 m_eigenvectorsOk(false),
126 m_realSchur(size),
127 m_matT(size, size),
128 m_tmp(size)
129 {}
130
131 /** \brief Constructor; computes eigendecomposition of given matrix.
132 *
133 * \param[in] matrix Square matrix whose eigendecomposition is to be computed.
134 * \param[in] computeEigenvectors If true, both the eigenvectors and the
135 * eigenvalues are computed; if false, only the eigenvalues are
136 * computed.
137 *
138 * This constructor calls compute() to compute the eigenvalues
139 * and eigenvectors.
140 *
141 * Example: \include EigenSolver_EigenSolver_MatrixType.cpp
142 * Output: \verbinclude EigenSolver_EigenSolver_MatrixType.out
143 *
144 * \sa compute()
145 */
146 template<typename InputType>
147 explicit EigenSolver(const EigenBase<InputType>& matrix, bool computeEigenvectors = true)
148 : m_eivec(matrix.rows(), matrix.cols()),
149 m_eivalues(matrix.cols()),
150 m_isInitialized(false),
151 m_eigenvectorsOk(false),
152 m_realSchur(matrix.cols()),
153 m_matT(matrix.rows(), matrix.cols()),
154 m_tmp(matrix.cols())
155 {
156 compute(matrix.derived(), computeEigenvectors);
157 }
158
159 /** \brief Returns the eigenvectors of given matrix.
160 *
161 * \returns %Matrix whose columns are the (possibly complex) eigenvectors.
162 *
163 * \pre Either the constructor
164 * EigenSolver(const MatrixType&,bool) or the member function
165 * compute(const MatrixType&, bool) has been called before, and
166 * \p computeEigenvectors was set to true (the default).
167 *
168 * Column \f$ k \f$ of the returned matrix is an eigenvector corresponding
169 * to eigenvalue number \f$ k \f$ as returned by eigenvalues(). The
170 * eigenvectors are normalized to have (Euclidean) norm equal to one. The
171 * matrix returned by this function is the matrix \f$ V \f$ in the
172 * eigendecomposition \f$ A = V D V^{-1} \f$, if it exists.
173 *
174 * Example: \include EigenSolver_eigenvectors.cpp
175 * Output: \verbinclude EigenSolver_eigenvectors.out
176 *
177 * \sa eigenvalues(), pseudoEigenvectors()
178 */
179 EigenvectorsType eigenvectors() const;
180
181 /** \brief Returns the pseudo-eigenvectors of given matrix.
182 *
183 * \returns Const reference to matrix whose columns are the pseudo-eigenvectors.
184 *
185 * \pre Either the constructor
186 * EigenSolver(const MatrixType&,bool) or the member function
187 * compute(const MatrixType&, bool) has been called before, and
188 * \p computeEigenvectors was set to true (the default).
189 *
190 * The real matrix \f$ V \f$ returned by this function and the
191 * block-diagonal matrix \f$ D \f$ returned by pseudoEigenvalueMatrix()
192 * satisfy \f$ AV = VD \f$.
193 *
194 * Example: \include EigenSolver_pseudoEigenvectors.cpp
195 * Output: \verbinclude EigenSolver_pseudoEigenvectors.out
196 *
197 * \sa pseudoEigenvalueMatrix(), eigenvectors()
198 */
pseudoEigenvectors()199 const MatrixType& pseudoEigenvectors() const
200 {
201 eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
202 eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
203 return m_eivec;
204 }
205
206 /** \brief Returns the block-diagonal matrix in the pseudo-eigendecomposition.
207 *
208 * \returns A block-diagonal matrix.
209 *
210 * \pre Either the constructor
211 * EigenSolver(const MatrixType&,bool) or the member function
212 * compute(const MatrixType&, bool) has been called before.
213 *
214 * The matrix \f$ D \f$ returned by this function is real and
215 * block-diagonal. The blocks on the diagonal are either 1-by-1 or 2-by-2
216 * blocks of the form
217 * \f$ \begin{bmatrix} u & v \\ -v & u \end{bmatrix} \f$.
218 * These blocks are not sorted in any particular order.
219 * The matrix \f$ D \f$ and the matrix \f$ V \f$ returned by
220 * pseudoEigenvectors() satisfy \f$ AV = VD \f$.
221 *
222 * \sa pseudoEigenvectors() for an example, eigenvalues()
223 */
224 MatrixType pseudoEigenvalueMatrix() const;
225
226 /** \brief Returns the eigenvalues of given matrix.
227 *
228 * \returns A const reference to the column vector containing the eigenvalues.
229 *
230 * \pre Either the constructor
231 * EigenSolver(const MatrixType&,bool) or the member function
232 * compute(const MatrixType&, bool) has been called before.
233 *
234 * The eigenvalues are repeated according to their algebraic multiplicity,
235 * so there are as many eigenvalues as rows in the matrix. The eigenvalues
236 * are not sorted in any particular order.
237 *
238 * Example: \include EigenSolver_eigenvalues.cpp
239 * Output: \verbinclude EigenSolver_eigenvalues.out
240 *
241 * \sa eigenvectors(), pseudoEigenvalueMatrix(),
242 * MatrixBase::eigenvalues()
243 */
eigenvalues()244 const EigenvalueType& eigenvalues() const
245 {
246 eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
247 return m_eivalues;
248 }
249
250 /** \brief Computes eigendecomposition of given matrix.
251 *
252 * \param[in] matrix Square matrix whose eigendecomposition is to be computed.
253 * \param[in] computeEigenvectors If true, both the eigenvectors and the
254 * eigenvalues are computed; if false, only the eigenvalues are
255 * computed.
256 * \returns Reference to \c *this
257 *
258 * This function computes the eigenvalues of the real matrix \p matrix.
259 * The eigenvalues() function can be used to retrieve them. If
260 * \p computeEigenvectors is true, then the eigenvectors are also computed
261 * and can be retrieved by calling eigenvectors().
262 *
263 * The matrix is first reduced to real Schur form using the RealSchur
264 * class. The Schur decomposition is then used to compute the eigenvalues
265 * and eigenvectors.
266 *
267 * The cost of the computation is dominated by the cost of the
268 * Schur decomposition, which is very approximately \f$ 25n^3 \f$
269 * (where \f$ n \f$ is the size of the matrix) if \p computeEigenvectors
270 * is true, and \f$ 10n^3 \f$ if \p computeEigenvectors is false.
271 *
272 * This method reuses of the allocated data in the EigenSolver object.
273 *
274 * Example: \include EigenSolver_compute.cpp
275 * Output: \verbinclude EigenSolver_compute.out
276 */
277 template<typename InputType>
278 EigenSolver& compute(const EigenBase<InputType>& matrix, bool computeEigenvectors = true);
279
280 /** \returns NumericalIssue if the input contains INF or NaN values or overflow occured. Returns Success otherwise. */
info()281 ComputationInfo info() const
282 {
283 eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
284 return m_info;
285 }
286
287 /** \brief Sets the maximum number of iterations allowed. */
setMaxIterations(Index maxIters)288 EigenSolver& setMaxIterations(Index maxIters)
289 {
290 m_realSchur.setMaxIterations(maxIters);
291 return *this;
292 }
293
294 /** \brief Returns the maximum number of iterations. */
getMaxIterations()295 Index getMaxIterations()
296 {
297 return m_realSchur.getMaxIterations();
298 }
299
300 private:
301 void doComputeEigenvectors();
302
303 protected:
304
check_template_parameters()305 static void check_template_parameters()
306 {
307 EIGEN_STATIC_ASSERT_NON_INTEGER(Scalar);
308 EIGEN_STATIC_ASSERT(!NumTraits<Scalar>::IsComplex, NUMERIC_TYPE_MUST_BE_REAL);
309 }
310
311 MatrixType m_eivec;
312 EigenvalueType m_eivalues;
313 bool m_isInitialized;
314 bool m_eigenvectorsOk;
315 ComputationInfo m_info;
316 RealSchur<MatrixType> m_realSchur;
317 MatrixType m_matT;
318
319 typedef Matrix<Scalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> ColumnVectorType;
320 ColumnVectorType m_tmp;
321 };
322
323 template<typename MatrixType>
pseudoEigenvalueMatrix()324 MatrixType EigenSolver<MatrixType>::pseudoEigenvalueMatrix() const
325 {
326 eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
327 const RealScalar precision = RealScalar(2)*NumTraits<RealScalar>::epsilon();
328 Index n = m_eivalues.rows();
329 MatrixType matD = MatrixType::Zero(n,n);
330 for (Index i=0; i<n; ++i)
331 {
332 if (internal::isMuchSmallerThan(numext::imag(m_eivalues.coeff(i)), numext::real(m_eivalues.coeff(i)), precision))
333 matD.coeffRef(i,i) = numext::real(m_eivalues.coeff(i));
334 else
335 {
336 matD.template block<2,2>(i,i) << numext::real(m_eivalues.coeff(i)), numext::imag(m_eivalues.coeff(i)),
337 -numext::imag(m_eivalues.coeff(i)), numext::real(m_eivalues.coeff(i));
338 ++i;
339 }
340 }
341 return matD;
342 }
343
344 template<typename MatrixType>
eigenvectors()345 typename EigenSolver<MatrixType>::EigenvectorsType EigenSolver<MatrixType>::eigenvectors() const
346 {
347 eigen_assert(m_isInitialized && "EigenSolver is not initialized.");
348 eigen_assert(m_eigenvectorsOk && "The eigenvectors have not been computed together with the eigenvalues.");
349 const RealScalar precision = RealScalar(2)*NumTraits<RealScalar>::epsilon();
350 Index n = m_eivec.cols();
351 EigenvectorsType matV(n,n);
352 for (Index j=0; j<n; ++j)
353 {
354 if (internal::isMuchSmallerThan(numext::imag(m_eivalues.coeff(j)), numext::real(m_eivalues.coeff(j)), precision) || j+1==n)
355 {
356 // we have a real eigen value
357 matV.col(j) = m_eivec.col(j).template cast<ComplexScalar>();
358 matV.col(j).normalize();
359 }
360 else
361 {
362 // we have a pair of complex eigen values
363 for (Index i=0; i<n; ++i)
364 {
365 matV.coeffRef(i,j) = ComplexScalar(m_eivec.coeff(i,j), m_eivec.coeff(i,j+1));
366 matV.coeffRef(i,j+1) = ComplexScalar(m_eivec.coeff(i,j), -m_eivec.coeff(i,j+1));
367 }
368 matV.col(j).normalize();
369 matV.col(j+1).normalize();
370 ++j;
371 }
372 }
373 return matV;
374 }
375
376 template<typename MatrixType>
377 template<typename InputType>
378 EigenSolver<MatrixType>&
compute(const EigenBase<InputType> & matrix,bool computeEigenvectors)379 EigenSolver<MatrixType>::compute(const EigenBase<InputType>& matrix, bool computeEigenvectors)
380 {
381 check_template_parameters();
382
383 using std::sqrt;
384 using std::abs;
385 using numext::isfinite;
386 eigen_assert(matrix.cols() == matrix.rows());
387
388 // Reduce to real Schur form.
389 m_realSchur.compute(matrix.derived(), computeEigenvectors);
390
391 m_info = m_realSchur.info();
392
393 if (m_info == Success)
394 {
395 m_matT = m_realSchur.matrixT();
396 if (computeEigenvectors)
397 m_eivec = m_realSchur.matrixU();
398
399 // Compute eigenvalues from matT
400 m_eivalues.resize(matrix.cols());
401 Index i = 0;
402 while (i < matrix.cols())
403 {
404 if (i == matrix.cols() - 1 || m_matT.coeff(i+1, i) == Scalar(0))
405 {
406 m_eivalues.coeffRef(i) = m_matT.coeff(i, i);
407 if(!(isfinite)(m_eivalues.coeffRef(i)))
408 {
409 m_isInitialized = true;
410 m_eigenvectorsOk = false;
411 m_info = NumericalIssue;
412 return *this;
413 }
414 ++i;
415 }
416 else
417 {
418 Scalar p = Scalar(0.5) * (m_matT.coeff(i, i) - m_matT.coeff(i+1, i+1));
419 Scalar z;
420 // Compute z = sqrt(abs(p * p + m_matT.coeff(i+1, i) * m_matT.coeff(i, i+1)));
421 // without overflow
422 {
423 Scalar t0 = m_matT.coeff(i+1, i);
424 Scalar t1 = m_matT.coeff(i, i+1);
425 Scalar maxval = numext::maxi<Scalar>(abs(p),numext::maxi<Scalar>(abs(t0),abs(t1)));
426 t0 /= maxval;
427 t1 /= maxval;
428 Scalar p0 = p/maxval;
429 z = maxval * sqrt(abs(p0 * p0 + t0 * t1));
430 }
431
432 m_eivalues.coeffRef(i) = ComplexScalar(m_matT.coeff(i+1, i+1) + p, z);
433 m_eivalues.coeffRef(i+1) = ComplexScalar(m_matT.coeff(i+1, i+1) + p, -z);
434 if(!((isfinite)(m_eivalues.coeffRef(i)) && (isfinite)(m_eivalues.coeffRef(i+1))))
435 {
436 m_isInitialized = true;
437 m_eigenvectorsOk = false;
438 m_info = NumericalIssue;
439 return *this;
440 }
441 i += 2;
442 }
443 }
444
445 // Compute eigenvectors.
446 if (computeEigenvectors)
447 doComputeEigenvectors();
448 }
449
450 m_isInitialized = true;
451 m_eigenvectorsOk = computeEigenvectors;
452
453 return *this;
454 }
455
456
457 template<typename MatrixType>
doComputeEigenvectors()458 void EigenSolver<MatrixType>::doComputeEigenvectors()
459 {
460 using std::abs;
461 const Index size = m_eivec.cols();
462 const Scalar eps = NumTraits<Scalar>::epsilon();
463
464 // inefficient! this is already computed in RealSchur
465 Scalar norm(0);
466 for (Index j = 0; j < size; ++j)
467 {
468 norm += m_matT.row(j).segment((std::max)(j-1,Index(0)), size-(std::max)(j-1,Index(0))).cwiseAbs().sum();
469 }
470
471 // Backsubstitute to find vectors of upper triangular form
472 if (norm == Scalar(0))
473 {
474 return;
475 }
476
477 for (Index n = size-1; n >= 0; n--)
478 {
479 Scalar p = m_eivalues.coeff(n).real();
480 Scalar q = m_eivalues.coeff(n).imag();
481
482 // Scalar vector
483 if (q == Scalar(0))
484 {
485 Scalar lastr(0), lastw(0);
486 Index l = n;
487
488 m_matT.coeffRef(n,n) = Scalar(1);
489 for (Index i = n-1; i >= 0; i--)
490 {
491 Scalar w = m_matT.coeff(i,i) - p;
492 Scalar r = m_matT.row(i).segment(l,n-l+1).dot(m_matT.col(n).segment(l, n-l+1));
493
494 if (m_eivalues.coeff(i).imag() < Scalar(0))
495 {
496 lastw = w;
497 lastr = r;
498 }
499 else
500 {
501 l = i;
502 if (m_eivalues.coeff(i).imag() == Scalar(0))
503 {
504 if (w != Scalar(0))
505 m_matT.coeffRef(i,n) = -r / w;
506 else
507 m_matT.coeffRef(i,n) = -r / (eps * norm);
508 }
509 else // Solve real equations
510 {
511 Scalar x = m_matT.coeff(i,i+1);
512 Scalar y = m_matT.coeff(i+1,i);
513 Scalar denom = (m_eivalues.coeff(i).real() - p) * (m_eivalues.coeff(i).real() - p) + m_eivalues.coeff(i).imag() * m_eivalues.coeff(i).imag();
514 Scalar t = (x * lastr - lastw * r) / denom;
515 m_matT.coeffRef(i,n) = t;
516 if (abs(x) > abs(lastw))
517 m_matT.coeffRef(i+1,n) = (-r - w * t) / x;
518 else
519 m_matT.coeffRef(i+1,n) = (-lastr - y * t) / lastw;
520 }
521
522 // Overflow control
523 Scalar t = abs(m_matT.coeff(i,n));
524 if ((eps * t) * t > Scalar(1))
525 m_matT.col(n).tail(size-i) /= t;
526 }
527 }
528 }
529 else if (q < Scalar(0) && n > 0) // Complex vector
530 {
531 Scalar lastra(0), lastsa(0), lastw(0);
532 Index l = n-1;
533
534 // Last vector component imaginary so matrix is triangular
535 if (abs(m_matT.coeff(n,n-1)) > abs(m_matT.coeff(n-1,n)))
536 {
537 m_matT.coeffRef(n-1,n-1) = q / m_matT.coeff(n,n-1);
538 m_matT.coeffRef(n-1,n) = -(m_matT.coeff(n,n) - p) / m_matT.coeff(n,n-1);
539 }
540 else
541 {
542 ComplexScalar cc = ComplexScalar(Scalar(0),-m_matT.coeff(n-1,n)) / ComplexScalar(m_matT.coeff(n-1,n-1)-p,q);
543 m_matT.coeffRef(n-1,n-1) = numext::real(cc);
544 m_matT.coeffRef(n-1,n) = numext::imag(cc);
545 }
546 m_matT.coeffRef(n,n-1) = Scalar(0);
547 m_matT.coeffRef(n,n) = Scalar(1);
548 for (Index i = n-2; i >= 0; i--)
549 {
550 Scalar ra = m_matT.row(i).segment(l, n-l+1).dot(m_matT.col(n-1).segment(l, n-l+1));
551 Scalar sa = m_matT.row(i).segment(l, n-l+1).dot(m_matT.col(n).segment(l, n-l+1));
552 Scalar w = m_matT.coeff(i,i) - p;
553
554 if (m_eivalues.coeff(i).imag() < Scalar(0))
555 {
556 lastw = w;
557 lastra = ra;
558 lastsa = sa;
559 }
560 else
561 {
562 l = i;
563 if (m_eivalues.coeff(i).imag() == RealScalar(0))
564 {
565 ComplexScalar cc = ComplexScalar(-ra,-sa) / ComplexScalar(w,q);
566 m_matT.coeffRef(i,n-1) = numext::real(cc);
567 m_matT.coeffRef(i,n) = numext::imag(cc);
568 }
569 else
570 {
571 // Solve complex equations
572 Scalar x = m_matT.coeff(i,i+1);
573 Scalar y = m_matT.coeff(i+1,i);
574 Scalar vr = (m_eivalues.coeff(i).real() - p) * (m_eivalues.coeff(i).real() - p) + m_eivalues.coeff(i).imag() * m_eivalues.coeff(i).imag() - q * q;
575 Scalar vi = (m_eivalues.coeff(i).real() - p) * Scalar(2) * q;
576 if ((vr == Scalar(0)) && (vi == Scalar(0)))
577 vr = eps * norm * (abs(w) + abs(q) + abs(x) + abs(y) + abs(lastw));
578
579 ComplexScalar cc = ComplexScalar(x*lastra-lastw*ra+q*sa,x*lastsa-lastw*sa-q*ra) / ComplexScalar(vr,vi);
580 m_matT.coeffRef(i,n-1) = numext::real(cc);
581 m_matT.coeffRef(i,n) = numext::imag(cc);
582 if (abs(x) > (abs(lastw) + abs(q)))
583 {
584 m_matT.coeffRef(i+1,n-1) = (-ra - w * m_matT.coeff(i,n-1) + q * m_matT.coeff(i,n)) / x;
585 m_matT.coeffRef(i+1,n) = (-sa - w * m_matT.coeff(i,n) - q * m_matT.coeff(i,n-1)) / x;
586 }
587 else
588 {
589 cc = ComplexScalar(-lastra-y*m_matT.coeff(i,n-1),-lastsa-y*m_matT.coeff(i,n)) / ComplexScalar(lastw,q);
590 m_matT.coeffRef(i+1,n-1) = numext::real(cc);
591 m_matT.coeffRef(i+1,n) = numext::imag(cc);
592 }
593 }
594
595 // Overflow control
596 Scalar t = numext::maxi<Scalar>(abs(m_matT.coeff(i,n-1)),abs(m_matT.coeff(i,n)));
597 if ((eps * t) * t > Scalar(1))
598 m_matT.block(i, n-1, size-i, 2) /= t;
599
600 }
601 }
602
603 // We handled a pair of complex conjugate eigenvalues, so need to skip them both
604 n--;
605 }
606 else
607 {
608 eigen_assert(0 && "Internal bug in EigenSolver (INF or NaN has not been detected)"); // this should not happen
609 }
610 }
611
612 // Back transformation to get eigenvectors of original matrix
613 for (Index j = size-1; j >= 0; j--)
614 {
615 m_tmp.noalias() = m_eivec.leftCols(j+1) * m_matT.col(j).segment(0, j+1);
616 m_eivec.col(j) = m_tmp;
617 }
618 }
619
620 } // end namespace Eigen
621
622 #endif // EIGEN_EIGENSOLVER_H
623