RealQZ.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2012 Alexey Korepanov <kaikaikai@yandex.ru>
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_REAL_QZ_H
11 #define EIGEN_REAL_QZ_H
12 
13 #include "./InternalHeaderCheck.h"
14 
15 namespace Eigen {
16 
59  template<typename MatrixType_> class RealQZ
60  {
61  public:
62  typedef MatrixType_ MatrixType;
63  enum {
64  RowsAtCompileTime = MatrixType::RowsAtCompileTime,
65  ColsAtCompileTime = MatrixType::ColsAtCompileTime,
66  Options = MatrixType::Options,
67  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
68  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
69  };
70  typedef typename MatrixType::Scalar Scalar;
71  typedef std::complex<typename NumTraits<Scalar>::Real> ComplexScalar;
72  typedef Eigen::Index Index;
73 
76 
89  m_S(size, size),
90  m_T(size, size),
91  m_Q(size, size),
92  m_Z(size, size),
93  m_workspace(size*2),
94  m_maxIters(400),
95  m_isInitialized(false),
96  m_computeQZ(true)
97  {}
98 
107  RealQZ(const MatrixType& A, const MatrixType& B, bool computeQZ = true) :
108  m_S(A.rows(),A.cols()),
109  m_T(A.rows(),A.cols()),
110  m_Q(A.rows(),A.cols()),
111  m_Z(A.rows(),A.cols()),
112  m_workspace(A.rows()*2),
113  m_maxIters(400),
114  m_isInitialized(false),
115  m_computeQZ(true)
116  {
117  compute(A, B, computeQZ);
118  }
119 
124  const MatrixType& matrixQ() const {
125  eigen_assert(m_isInitialized && "RealQZ is not initialized.");
126  eigen_assert(m_computeQZ && "The matrices Q and Z have not been computed during the QZ decomposition.");
127  return m_Q;
128  }
129 
134  const MatrixType& matrixZ() const {
135  eigen_assert(m_isInitialized && "RealQZ is not initialized.");
136  eigen_assert(m_computeQZ && "The matrices Q and Z have not been computed during the QZ decomposition.");
137  return m_Z;
138  }
139 
144  const MatrixType& matrixS() const {
145  eigen_assert(m_isInitialized && "RealQZ is not initialized.");
146  return m_S;
147  }
148 
153  const MatrixType& matrixT() const {
154  eigen_assert(m_isInitialized && "RealQZ is not initialized.");
155  return m_T;
156  }
157 
165  RealQZ& compute(const MatrixType& A, const MatrixType& B, bool computeQZ = true);
166 
172  {
173  eigen_assert(m_isInitialized && "RealQZ is not initialized.");
174  return m_info;
175  }
176 
180  {
181  eigen_assert(m_isInitialized && "RealQZ is not initialized.");
182  return m_global_iter;
183  }
184 
189  {
190  m_maxIters = maxIters;
191  return *this;
192  }
193 
194  private:
195 
204 
209 
211  void computeNorms();
215  void pushDownZero(Index z, Index f, Index l);
216  void step(Index f, Index l, Index iter);
217 
218  }; // RealQZ
219 
221  template<typename MatrixType>
223  {
224 
225  const Index dim = m_S.cols();
226 
227  // perform QR decomposition of T, overwrite T with R, save Q
228  HouseholderQR<MatrixType> qrT(m_T);
229  m_T = qrT.matrixQR();
230  m_T.template triangularView<StrictlyLower>().setZero();
231  m_Q = qrT.householderQ();
232  // overwrite S with Q* S
233  m_S.applyOnTheLeft(m_Q.adjoint());
234  // init Z as Identity
235  if (m_computeQZ)
236  m_Z = MatrixType::Identity(dim,dim);
237  // reduce S to upper Hessenberg with Givens rotations
238  for (Index j=0; j<=dim-3; j++) {
239  for (Index i=dim-1; i>=j+2; i--) {
240  JRs G;
241  // kill S(i,j)
242  if(!numext::is_exactly_zero(m_S.coeff(i, j)))
243  {
244  G.makeGivens(m_S.coeff(i-1,j), m_S.coeff(i,j), &m_S.coeffRef(i-1, j));
245  m_S.coeffRef(i,j) = Scalar(0.0);
246  m_S.rightCols(dim-j-1).applyOnTheLeft(i-1,i,G.adjoint());
247  m_T.rightCols(dim-i+1).applyOnTheLeft(i-1,i,G.adjoint());
248  // update Q
249  if (m_computeQZ)
250  m_Q.applyOnTheRight(i-1,i,G);
251  }
252  // kill T(i,i-1)
253  if(!numext::is_exactly_zero(m_T.coeff(i, i - 1)))
254  {
255  G.makeGivens(m_T.coeff(i,i), m_T.coeff(i,i-1), &m_T.coeffRef(i,i));
256  m_T.coeffRef(i,i-1) = Scalar(0.0);
257  m_S.applyOnTheRight(i,i-1,G);
258  m_T.topRows(i).applyOnTheRight(i,i-1,G);
259  // update Z
260  if (m_computeQZ)
261  m_Z.applyOnTheLeft(i,i-1,G.adjoint());
262  }
263  }
264  }
265  }
266 
268  template<typename MatrixType>
270  {
271  const Index size = m_S.cols();
272  m_normOfS = Scalar(0.0);
273  m_normOfT = Scalar(0.0);
274  for (Index j = 0; j < size; ++j)
275  {
276  m_normOfS += m_S.col(j).segment(0, (std::min)(size,j+2)).cwiseAbs().sum();
277  m_normOfT += m_T.row(j).segment(j, size - j).cwiseAbs().sum();
278  }
279  }
280 
281 
283  template<typename MatrixType>
285  {
286  using std::abs;
287  Index res = iu;
288  while (res > 0)
289  {
290  Scalar s = abs(m_S.coeff(res-1,res-1)) + abs(m_S.coeff(res,res));
292  s = m_normOfS;
293  if (abs(m_S.coeff(res,res-1)) < NumTraits<Scalar>::epsilon() * s)
294  break;
295  res--;
296  }
297  return res;
298  }
299 
301  template<typename MatrixType>
303  {
304  using std::abs;
305  Index res = l;
306  while (res >= f) {
307  if (abs(m_T.coeff(res,res)) <= NumTraits<Scalar>::epsilon() * m_normOfT)
308  break;
309  res--;
310  }
311  return res;
312  }
313 
315  template<typename MatrixType>
317  {
318  using std::abs;
319  using std::sqrt;
320  const Index dim=m_S.cols();
321  if (numext::is_exactly_zero(abs(m_S.coeff(i + 1, i))))
322  return;
323  Index j = findSmallDiagEntry(i,i+1);
324  if (j==i-1)
325  {
326  // block of (S T^{-1})
327  Matrix2s STi = m_T.template block<2,2>(i,i).template triangularView<Upper>().
328  template solve<OnTheRight>(m_S.template block<2,2>(i,i));
329  Scalar p = Scalar(0.5)*(STi(0,0)-STi(1,1));
330  Scalar q = p*p + STi(1,0)*STi(0,1);
331  if (q>=0) {
332  Scalar z = sqrt(q);
333  // one QR-like iteration for ABi - lambda I
334  // is enough - when we know exact eigenvalue in advance,
335  // convergence is immediate
336  JRs G;
337  if (p>=0)
338  G.makeGivens(p + z, STi(1,0));
339  else
340  G.makeGivens(p - z, STi(1,0));
341  m_S.rightCols(dim-i).applyOnTheLeft(i,i+1,G.adjoint());
342  m_T.rightCols(dim-i).applyOnTheLeft(i,i+1,G.adjoint());
343  // update Q
344  if (m_computeQZ)
345  m_Q.applyOnTheRight(i,i+1,G);
346 
347  G.makeGivens(m_T.coeff(i+1,i+1), m_T.coeff(i+1,i));
348  m_S.topRows(i+2).applyOnTheRight(i+1,i,G);
349  m_T.topRows(i+2).applyOnTheRight(i+1,i,G);
350  // update Z
351  if (m_computeQZ)
352  m_Z.applyOnTheLeft(i+1,i,G.adjoint());
353 
354  m_S.coeffRef(i+1,i) = Scalar(0.0);
355  m_T.coeffRef(i+1,i) = Scalar(0.0);
356  }
357  }
358  else
359  {
360  pushDownZero(j,i,i+1);
361  }
362  }
363 
365  template<typename MatrixType>
367  {
368  JRs G;
369  const Index dim = m_S.cols();
370  for (Index zz=z; zz<l; zz++)
371  {
372  // push 0 down
373  Index firstColS = zz>f ? (zz-1) : zz;
374  G.makeGivens(m_T.coeff(zz, zz+1), m_T.coeff(zz+1, zz+1));
375  m_S.rightCols(dim-firstColS).applyOnTheLeft(zz,zz+1,G.adjoint());
376  m_T.rightCols(dim-zz).applyOnTheLeft(zz,zz+1,G.adjoint());
377  m_T.coeffRef(zz+1,zz+1) = Scalar(0.0);
378  // update Q
379  if (m_computeQZ)
380  m_Q.applyOnTheRight(zz,zz+1,G);
381  // kill S(zz+1, zz-1)
382  if (zz>f)
383  {
384  G.makeGivens(m_S.coeff(zz+1, zz), m_S.coeff(zz+1,zz-1));
385  m_S.topRows(zz+2).applyOnTheRight(zz, zz-1,G);
386  m_T.topRows(zz+1).applyOnTheRight(zz, zz-1,G);
387  m_S.coeffRef(zz+1,zz-1) = Scalar(0.0);
388  // update Z
389  if (m_computeQZ)
390  m_Z.applyOnTheLeft(zz,zz-1,G.adjoint());
391  }
392  }
393  // finally kill S(l,l-1)
394  G.makeGivens(m_S.coeff(l,l), m_S.coeff(l,l-1));
395  m_S.applyOnTheRight(l,l-1,G);
396  m_T.applyOnTheRight(l,l-1,G);
397  m_S.coeffRef(l,l-1)=Scalar(0.0);
398  // update Z
399  if (m_computeQZ)
400  m_Z.applyOnTheLeft(l,l-1,G.adjoint());
401  }
402 
404  template<typename MatrixType>
405  inline void RealQZ<MatrixType>::step(Index f, Index l, Index iter)
406  {
407  using std::abs;
408  const Index dim = m_S.cols();
409 
410  // x, y, z
411  Scalar x, y, z;
412  if (iter==10)
413  {
414  // Wilkinson ad hoc shift
415  const Scalar
416  a11=m_S.coeff(f+0,f+0), a12=m_S.coeff(f+0,f+1),
417  a21=m_S.coeff(f+1,f+0), a22=m_S.coeff(f+1,f+1), a32=m_S.coeff(f+2,f+1),
418  b12=m_T.coeff(f+0,f+1),
419  b11i=Scalar(1.0)/m_T.coeff(f+0,f+0),
420  b22i=Scalar(1.0)/m_T.coeff(f+1,f+1),
421  a87=m_S.coeff(l-1,l-2),
422  a98=m_S.coeff(l-0,l-1),
423  b77i=Scalar(1.0)/m_T.coeff(l-2,l-2),
424  b88i=Scalar(1.0)/m_T.coeff(l-1,l-1);
425  Scalar ss = abs(a87*b77i) + abs(a98*b88i),
426  lpl = Scalar(1.5)*ss,
427  ll = ss*ss;
428  x = ll + a11*a11*b11i*b11i - lpl*a11*b11i + a12*a21*b11i*b22i
429  - a11*a21*b12*b11i*b11i*b22i;
430  y = a11*a21*b11i*b11i - lpl*a21*b11i + a21*a22*b11i*b22i
431  - a21*a21*b12*b11i*b11i*b22i;
432  z = a21*a32*b11i*b22i;
433  }
434  else if (iter==16)
435  {
436  // another exceptional shift
437  x = m_S.coeff(f,f)/m_T.coeff(f,f)-m_S.coeff(l,l)/m_T.coeff(l,l) + m_S.coeff(l,l-1)*m_T.coeff(l-1,l) /
438  (m_T.coeff(l-1,l-1)*m_T.coeff(l,l));
439  y = m_S.coeff(f+1,f)/m_T.coeff(f,f);
440  z = 0;
441  }
442  else if (iter>23 && !(iter%8))
443  {
444  // extremely exceptional shift
445  x = internal::random<Scalar>(-1.0,1.0);
446  y = internal::random<Scalar>(-1.0,1.0);
447  z = internal::random<Scalar>(-1.0,1.0);
448  }
449  else
450  {
451  // Compute the shifts: (x,y,z,0...) = (AB^-1 - l1 I) (AB^-1 - l2 I) e1
452  // where l1 and l2 are the eigenvalues of the 2x2 matrix C = U V^-1 where
453  // U and V are 2x2 bottom right sub matrices of A and B. Thus:
454  // = AB^-1AB^-1 + l1 l2 I - (l1+l2)(AB^-1)
455  // = AB^-1AB^-1 + det(M) - tr(M)(AB^-1)
456  // Since we are only interested in having x, y, z with a correct ratio, we have:
457  const Scalar
458  a11 = m_S.coeff(f,f), a12 = m_S.coeff(f,f+1),
459  a21 = m_S.coeff(f+1,f), a22 = m_S.coeff(f+1,f+1),
460  a32 = m_S.coeff(f+2,f+1),
461 
462  a88 = m_S.coeff(l-1,l-1), a89 = m_S.coeff(l-1,l),
463  a98 = m_S.coeff(l,l-1), a99 = m_S.coeff(l,l),
464 
465  b11 = m_T.coeff(f,f), b12 = m_T.coeff(f,f+1),
466  b22 = m_T.coeff(f+1,f+1),
467 
468  b88 = m_T.coeff(l-1,l-1), b89 = m_T.coeff(l-1,l),
469  b99 = m_T.coeff(l,l);
470 
471  x = ( (a88/b88 - a11/b11)*(a99/b99 - a11/b11) - (a89/b99)*(a98/b88) + (a98/b88)*(b89/b99)*(a11/b11) ) * (b11/a21)
472  + a12/b22 - (a11/b11)*(b12/b22);
473  y = (a22/b22-a11/b11) - (a21/b11)*(b12/b22) - (a88/b88-a11/b11) - (a99/b99-a11/b11) + (a98/b88)*(b89/b99);
474  z = a32/b22;
475  }
476 
477  JRs G;
478 
479  for (Index k=f; k<=l-2; k++)
480  {
481  // variables for Householder reflections
482  Vector2s essential2;
483  Scalar tau, beta;
484 
485  Vector3s hr(x,y,z);
486 
487  // Q_k to annihilate S(k+1,k-1) and S(k+2,k-1)
488  hr.makeHouseholderInPlace(tau, beta);
489  essential2 = hr.template bottomRows<2>();
490  Index fc=(std::max)(k-1,Index(0)); // first col to update
491  m_S.template middleRows<3>(k).rightCols(dim-fc).applyHouseholderOnTheLeft(essential2, tau, m_workspace.data());
492  m_T.template middleRows<3>(k).rightCols(dim-fc).applyHouseholderOnTheLeft(essential2, tau, m_workspace.data());
493  if (m_computeQZ)
494  m_Q.template middleCols<3>(k).applyHouseholderOnTheRight(essential2, tau, m_workspace.data());
495  if (k>f)
496  m_S.coeffRef(k+2,k-1) = m_S.coeffRef(k+1,k-1) = Scalar(0.0);
497 
498  // Z_{k1} to annihilate T(k+2,k+1) and T(k+2,k)
499  hr << m_T.coeff(k+2,k+2),m_T.coeff(k+2,k),m_T.coeff(k+2,k+1);
500  hr.makeHouseholderInPlace(tau, beta);
501  essential2 = hr.template bottomRows<2>();
502  {
503  Index lr = (std::min)(k+4,dim); // last row to update
504  Map<Matrix<Scalar,Dynamic,1> > tmp(m_workspace.data(),lr);
505  // S
506  tmp = m_S.template middleCols<2>(k).topRows(lr) * essential2;
507  tmp += m_S.col(k+2).head(lr);
508  m_S.col(k+2).head(lr) -= tau*tmp;
509  m_S.template middleCols<2>(k).topRows(lr) -= (tau*tmp) * essential2.adjoint();
510  // T
511  tmp = m_T.template middleCols<2>(k).topRows(lr) * essential2;
512  tmp += m_T.col(k+2).head(lr);
513  m_T.col(k+2).head(lr) -= tau*tmp;
514  m_T.template middleCols<2>(k).topRows(lr) -= (tau*tmp) * essential2.adjoint();
515  }
516  if (m_computeQZ)
517  {
518  // Z
519  Map<Matrix<Scalar,1,Dynamic> > tmp(m_workspace.data(),dim);
520  tmp = essential2.adjoint()*(m_Z.template middleRows<2>(k));
521  tmp += m_Z.row(k+2);
522  m_Z.row(k+2) -= tau*tmp;
523  m_Z.template middleRows<2>(k) -= essential2 * (tau*tmp);
524  }
525  m_T.coeffRef(k+2,k) = m_T.coeffRef(k+2,k+1) = Scalar(0.0);
526 
527  // Z_{k2} to annihilate T(k+1,k)
528  G.makeGivens(m_T.coeff(k+1,k+1), m_T.coeff(k+1,k));
529  m_S.applyOnTheRight(k+1,k,G);
530  m_T.applyOnTheRight(k+1,k,G);
531  // update Z
532  if (m_computeQZ)
533  m_Z.applyOnTheLeft(k+1,k,G.adjoint());
534  m_T.coeffRef(k+1,k) = Scalar(0.0);
535 
536  // update x,y,z
537  x = m_S.coeff(k+1,k);
538  y = m_S.coeff(k+2,k);
539  if (k < l-2)
540  z = m_S.coeff(k+3,k);
541  } // loop over k
542 
543  // Q_{n-1} to annihilate y = S(l,l-2)
544  G.makeGivens(x,y);
545  m_S.applyOnTheLeft(l-1,l,G.adjoint());
546  m_T.applyOnTheLeft(l-1,l,G.adjoint());
547  if (m_computeQZ)
548  m_Q.applyOnTheRight(l-1,l,G);
549  m_S.coeffRef(l,l-2) = Scalar(0.0);
550 
551  // Z_{n-1} to annihilate T(l,l-1)
552  G.makeGivens(m_T.coeff(l,l),m_T.coeff(l,l-1));
553  m_S.applyOnTheRight(l,l-1,G);
554  m_T.applyOnTheRight(l,l-1,G);
555  if (m_computeQZ)
556  m_Z.applyOnTheLeft(l,l-1,G.adjoint());
557  m_T.coeffRef(l,l-1) = Scalar(0.0);
558  }
559 
560  template<typename MatrixType>
561  RealQZ<MatrixType>& RealQZ<MatrixType>::compute(const MatrixType& A_in, const MatrixType& B_in, bool computeQZ)
562  {
563 
564  const Index dim = A_in.cols();
565 
566  eigen_assert (A_in.rows()==dim && A_in.cols()==dim
567  && B_in.rows()==dim && B_in.cols()==dim
568  && "Need square matrices of the same dimension");
569 
570  m_isInitialized = true;
571  m_computeQZ = computeQZ;
572  m_S = A_in; m_T = B_in;
573  m_workspace.resize(dim*2);
574  m_global_iter = 0;
575 
576  // entrance point: hessenberg triangular decomposition
577  hessenbergTriangular();
578  // compute L1 vector norms of T, S into m_normOfS, m_normOfT
579  computeNorms();
580 
581  Index l = dim-1,
582  f,
583  local_iter = 0;
584 
585  while (l>0 && local_iter<m_maxIters)
586  {
587  f = findSmallSubdiagEntry(l);
588  // now rows and columns f..l (including) decouple from the rest of the problem
589  if (f>0) m_S.coeffRef(f,f-1) = Scalar(0.0);
590  if (f == l) // One root found
591  {
592  l--;
593  local_iter = 0;
594  }
595  else if (f == l-1) // Two roots found
596  {
597  splitOffTwoRows(f);
598  l -= 2;
599  local_iter = 0;
600  }
601  else // No convergence yet
602  {
603  // if there's zero on diagonal of T, we can isolate an eigenvalue with Givens rotations
604  Index z = findSmallDiagEntry(f,l);
605  if (z>=f)
606  {
607  // zero found
608  pushDownZero(z,f,l);
609  }
610  else
611  {
612  // We are sure now that S.block(f,f, l-f+1,l-f+1) is underuced upper-Hessenberg
613  // and T.block(f,f, l-f+1,l-f+1) is invertible uper-triangular, which allows to
614  // apply a QR-like iteration to rows and columns f..l.
615  step(f,l, local_iter);
616  local_iter++;
617  m_global_iter++;
618  }
619  }
620  }
621  // check if we converged before reaching iterations limit
622  m_info = (local_iter<m_maxIters) ? Success : NoConvergence;
623 
624  // For each non triangular 2x2 diagonal block of S,
625  // reduce the respective 2x2 diagonal block of T to positive diagonal form using 2x2 SVD.
626  // This step is not mandatory for QZ, but it does help further extraction of eigenvalues/eigenvectors,
627  // and is in par with Lapack/Matlab QZ.
628  if(m_info==Success)
629  {
630  for(Index i=0; i<dim-1; ++i)
631  {
632  if(!numext::is_exactly_zero(m_S.coeff(i + 1, i)))
633  {
634  JacobiRotation<Scalar> j_left, j_right;
635  internal::real_2x2_jacobi_svd(m_T, i, i+1, &j_left, &j_right);
636 
637  // Apply resulting Jacobi rotations
638  m_S.applyOnTheLeft(i,i+1,j_left);
639  m_S.applyOnTheRight(i,i+1,j_right);
640  m_T.applyOnTheLeft(i,i+1,j_left);
641  m_T.applyOnTheRight(i,i+1,j_right);
642  m_T(i+1,i) = m_T(i,i+1) = Scalar(0);
643 
644  if(m_computeQZ) {
645  m_Q.applyOnTheRight(i,i+1,j_left.transpose());
646  m_Z.applyOnTheLeft(i,i+1,j_right.transpose());
647  }
648 
649  i++;
650  }
651  }
652  }
653 
654  return *this;
655  } // end compute
656 
657 } // end namespace Eigen
658 
659 #endif //EIGEN_REAL_QZ
const AbsReturnType abs() const
const SqrtReturnType sqrt() const
solver compute(A)
MatrixXcf A
MatrixXf B
JacobiRotation< float > G
#define eigen_assert(x)
Definition: Macros.h:902
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
float * p
internal::traits< Derived >::Scalar Scalar
Definition: DenseBase.h:61
Householder QR decomposition of a matrix.
Definition: HouseholderQR.h:60
HouseholderSequenceType householderQ() const
const MatrixType & matrixQR() const
Rotation given by a cosine-sine pair.
Definition: Jacobi.h:37
JacobiRotation transpose() const
Definition: Jacobi.h:65
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:98
const AdjointReturnType adjoint() const
Definition: Transpose.h:223
constexpr const Scalar & coeff(Index rowId, Index colId) const
const Scalar * data() const
constexpr Scalar & coeffRef(Index rowId, Index colId)
Performs a real QZ decomposition of a pair of square matrices.
Definition: RealQZ.h:60
void computeNorms()
Definition: RealQZ.h:269
@ ColsAtCompileTime
Definition: RealQZ.h:65
@ MaxRowsAtCompileTime
Definition: RealQZ.h:67
@ RowsAtCompileTime
Definition: RealQZ.h:64
@ MaxColsAtCompileTime
Definition: RealQZ.h:68
RealQZ & compute(const MatrixType &A, const MatrixType &B, bool computeQZ=true)
Computes QZ decomposition of given matrix.
Definition: RealQZ.h:561
Index iterations() const
Returns number of performed QR-like iterations.
Definition: RealQZ.h:179
Index m_global_iter
Definition: RealQZ.h:203
Matrix< Scalar, Dynamic, 1 > m_workspace
Definition: RealQZ.h:197
MatrixType_ MatrixType
Definition: RealQZ.h:62
Index findSmallSubdiagEntry(Index iu)
Definition: RealQZ.h:284
const MatrixType & matrixQ() const
Returns matrix Q in the QZ decomposition.
Definition: RealQZ.h:124
MatrixType m_Q
Definition: RealQZ.h:196
bool m_computeQZ
Definition: RealQZ.h:201
Matrix< Scalar, 2, 1 > Vector2s
Definition: RealQZ.h:206
RealQZ & setMaxIterations(Index maxIters)
Definition: RealQZ.h:188
RealQZ(const MatrixType &A, const MatrixType &B, bool computeQZ=true)
Constructor; computes real QZ decomposition of given matrices.
Definition: RealQZ.h:107
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: RealQZ.h:171
void splitOffTwoRows(Index i)
Definition: RealQZ.h:316
Scalar m_normOfS
Definition: RealQZ.h:202
const MatrixType & matrixZ() const
Returns matrix Z in the QZ decomposition.
Definition: RealQZ.h:134
Matrix< Scalar, 3, 1 > Vector3s
Definition: RealQZ.h:205
const MatrixType & matrixT() const
Returns matrix S in the QZ decomposition.
Definition: RealQZ.h:153
MatrixType::Scalar Scalar
Definition: RealQZ.h:70
void step(Index f, Index l, Index iter)
Definition: RealQZ.h:405
Matrix< Scalar, ColsAtCompileTime, 1, Options &~RowMajor, MaxColsAtCompileTime, 1 > ColumnVectorType
Definition: RealQZ.h:75
RealQZ(Index size=RowsAtCompileTime==Dynamic ? 1 :RowsAtCompileTime)
Default constructor.
Definition: RealQZ.h:88
void hessenbergTriangular()
Definition: RealQZ.h:222
MatrixType m_S
Definition: RealQZ.h:196
Matrix< ComplexScalar, ColsAtCompileTime, 1, Options &~RowMajor, MaxColsAtCompileTime, 1 > EigenvalueType
Definition: RealQZ.h:74
const MatrixType & matrixS() const
Returns matrix S in the QZ decomposition.
Definition: RealQZ.h:144
JacobiRotation< Scalar > JRs
Definition: RealQZ.h:208
Eigen::Index Index
Definition: RealQZ.h:72
Matrix< Scalar, 2, 2 > Matrix2s
Definition: RealQZ.h:207
ComputationInfo m_info
Definition: RealQZ.h:198
std::complex< typename NumTraits< Scalar >::Real > ComplexScalar
Definition: RealQZ.h:71
MatrixType m_Z
Definition: RealQZ.h:196
MatrixType m_T
Definition: RealQZ.h:196
bool m_isInitialized
Definition: RealQZ.h:200
Index m_maxIters
Definition: RealQZ.h:199
Index findSmallDiagEntry(Index f, Index l)
Definition: RealQZ.h:302
void pushDownZero(Index z, Index f, Index l)
Definition: RealQZ.h:366
Scalar m_normOfT
Definition: RealQZ.h:202
ComputationInfo
Definition: Constants.h:444
@ Success
Definition: Constants.h:446
@ NoConvergence
Definition: Constants.h:450
bfloat16() max(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:690
bfloat16() min(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:684
const Scalar & y
void real_2x2_jacobi_svd(const MatrixType &matrix, Index p, Index q, JacobiRotation< RealScalar > *j_left, JacobiRotation< RealScalar > *j_right)
Definition: RealSvd2x2.h:21
bool is_exactly_zero(const X &x)
Definition: Meta.h:475
: InteropHeaders
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
const int Dynamic
Definition: Constants.h:24
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:231
std::ptrdiff_t j