14#ifndef EIGEN_COMPLEX_QZ_H_
15#define EIGEN_COMPLEX_QZ_H_
18#include "./InternalHeaderCheck.h"
50template <
typename MatrixType_>
53 using MatrixType = MatrixType_;
54 using Scalar =
typename MatrixType_::Scalar;
55 using RealScalar =
typename MatrixType_::RealScalar;
58 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
59 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
60 Options = internal::traits<MatrixType>::Options,
61 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
62 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
65 using Vec = Matrix<Scalar, Dynamic, 1>;
66 using Vec2 = Matrix<Scalar, 2, 1>;
67 using Vec3 = Matrix<Scalar, 3, 1>;
68 using Row2 = Matrix<Scalar, 1, 2>;
69 using Mat2 = Matrix<Scalar, 2, 2>;
75 const MatrixType& matrixQ()
const {
76 eigen_assert(m_isInitialized &&
"ComplexQZ is not initialized.");
77 eigen_assert(m_computeQZ &&
"The matrices Q and Z have not been computed during the QZ decomposition.");
85 const MatrixType& matrixZ()
const {
86 eigen_assert(m_isInitialized &&
"ComplexQZ is not initialized.");
87 eigen_assert(m_computeQZ &&
"The matrices Q and Z have not been computed during the QZ decomposition.");
95 const MatrixType& matrixS()
const {
96 eigen_assert(m_isInitialized &&
"ComplexQZ is not initialized.");
104 const MatrixType& matrixT()
const {
105 eigen_assert(m_isInitialized &&
"ComplexQZ is not initialized.");
117 ComplexQZ(Index n,
bool computeQZ =
true,
unsigned int maxIters = 400)
121 m_Q(computeQZ ? n : (MatrixType::RowsAtCompileTime ==
Eigen::Dynamic ? 0 : MatrixType::RowsAtCompileTime),
122 computeQZ ? n : (MatrixType::ColsAtCompileTime ==
Eigen::Dynamic ? 0 : MatrixType::ColsAtCompileTime)),
123 m_Z(computeQZ ? n : (MatrixType::RowsAtCompileTime ==
Eigen::Dynamic ? 0 : MatrixType::RowsAtCompileTime),
124 computeQZ ? n : (MatrixType::ColsAtCompileTime ==
Eigen::Dynamic ? 0 : MatrixType::ColsAtCompileTime)),
126 m_computeQZ(computeQZ),
127 m_maxIters(maxIters) {}
140 ComplexQZ(
const MatrixType& A,
const MatrixType& B,
bool computeQZ =
true,
unsigned int maxIters = 400)
142 m_maxIters(maxIters),
143 m_computeQZ(computeQZ),
144 m_S(A.rows(), A.cols()),
145 m_T(A.rows(), A.cols()),
146 m_Q(computeQZ ? m_n : (MatrixType::RowsAtCompileTime ==
Eigen::Dynamic ? 0 : MatrixType::RowsAtCompileTime),
147 computeQZ ? m_n : (MatrixType::ColsAtCompileTime ==
Eigen::Dynamic ? 0 : MatrixType::ColsAtCompileTime)),
148 m_Z(computeQZ ? m_n : (MatrixType::RowsAtCompileTime ==
Eigen::Dynamic ? 0 : MatrixType::RowsAtCompileTime),
149 computeQZ ? m_n : (MatrixType::ColsAtCompileTime ==
Eigen::Dynamic ? 0 : MatrixType::ColsAtCompileTime)),
151 compute(A, B, computeQZ);
160 void compute(
const MatrixType& A,
const MatrixType& B,
bool computeQZ =
true);
170 template <
typename SparseMatrixType_>
171 void computeSparse(
const SparseMatrixType_& A,
const SparseMatrixType_& B,
bool computeQZ =
true);
177 ComputationInfo info()
const {
return m_info; }
181 unsigned int iterations()
const {
182 eigen_assert(m_isInitialized &&
"ComplexQZ is not initialized.");
183 return m_global_iter;
188 const unsigned int m_maxIters;
189 unsigned int m_global_iter;
190 bool m_isInitialized;
192 ComputationInfo m_info;
193 MatrixType m_S, m_T, m_Q, m_Z;
194 RealScalar m_normOfT, m_normOfS;
198 static bool is_negligible(
const Scalar x,
const RealScalar tol = NumTraits<RealScalar>::epsilon()) {
199 return numext::abs(x) <= tol;
202 void do_QZ_step(Index p, Index q);
204 inline Mat2 computeZk2(
const Row2& b);
207 void hessenbergTriangular(
const MatrixType& A,
const MatrixType& B);
211 void reduceHessenbergTriangular();
214 template <
typename SparseMatrixType_>
215 void hessenbergTriangularSparse(
const SparseMatrixType_& A,
const SparseMatrixType_& B);
219 Index findSmallSubdiagEntry(Index l);
220 Index findSmallDiagEntry(Index f, Index l);
222 void push_down_zero_ST(Index k, Index l);
224 void reduceDiagonal2x2block(Index i);
227template <
typename MatrixType_>
229 m_computeQZ = computeQZ;
232 eigen_assert(m_n == A.cols() &&
"A is not a square matrix");
233 eigen_assert(m_n == B.rows() && m_n == B.cols() &&
"B is not a square matrix or B is not of the same size as A");
235 m_isInitialized =
true;
239 hessenbergTriangular(A, B);
243 reduceHessenbergTriangular();
247template <
typename MatrixType_>
254 HouseholderQR<MatrixType> qr(m_T);
256 m_T.template triangularView<StrictlyLower>().setZero();
258 if (m_computeQZ) m_Q = qr.householderQ();
261 m_S.applyOnTheLeft(qr.householderQ().adjoint());
263 if (m_computeQZ) m_Z = MatrixType::Identity(m_n, m_n);
266 for (Index j = 0; j <= m_n - 3; j++) {
267 for (Index i = m_n - 1; i >= j + 2; i--) {
268 JacobiRotation<Scalar> G;
270 if (!numext::is_exactly_zero(m_S.coeff(i, j))) {
271 G.makeGivens(m_S.coeff(i - 1, j), m_S.coeff(i, j), &m_S.coeffRef(i - 1, j));
272 m_S.coeffRef(i, j) = Scalar(0);
273 m_T.rightCols(m_n - i + 1).applyOnTheLeft(i - 1, i, G.adjoint());
274 m_S.rightCols(m_n - j - 1).applyOnTheLeft(i - 1, i, G.adjoint());
276 if (!is_negligible(m_S(i, j)))
279 m_S(i, j) = Scalar(0);
281 if (m_computeQZ) m_Q.applyOnTheRight(i - 1, i, G);
284 if (!numext::is_exactly_zero(m_T.coeff(i, i - 1))) {
286 G.makeGivens(m_T.coeff(i, i), m_T.coeff(i, i - 1), &m_T.coeffRef(i, i));
287 m_T.topRows(i).applyOnTheRight(i - 1, i, G.adjoint());
288 m_T.coeffRef(i, i - 1) = Scalar(0);
290 m_S.applyOnTheRight(i - 1, i, G.adjoint());
292 if (m_computeQZ) m_Z.applyOnTheLeft(i - 1, i, G);
298template <
typename MatrixType>
299template <
typename SparseMatrixType_>
303 SparseQR<SparseMatrix<Scalar, ColMajor>, NaturalOrdering<Index>> sparseQR;
305 eigen_assert(B.isCompressed() &&
306 "SparseQR requires a sparse matrix in compressed mode."
307 "Call .makeCompressed() before passing it to SparseQR");
310 sparseQR.setPivotThreshold(RealScalar(0));
314 m_T = sparseQR.matrixR();
315 m_T.template triangularView<StrictlyLower>().setZero();
317 if (m_computeQZ) m_Q = sparseQR.matrixQ();
320 m_S = sparseQR.matrixQ().adjoint() * m_S;
322 if (m_computeQZ) m_Z = MatrixType::Identity(m_n, m_n);
325 for (Index j = 0; j <= m_n - 3; j++) {
326 for (Index i = m_n - 1; i >= j + 2; i--) {
327 JacobiRotation<Scalar> G;
330 if (m_S.coeff(i, j) != Scalar(0)) {
332 G.makeGivens(m_S.coeff(i - 1, j), m_S.coeff(i, j), &m_S.coeffRef(i - 1, j));
333 m_S.coeffRef(i, j) = Scalar(0);
334 m_T.rightCols(m_n - i + 1).applyOnTheLeft(i - 1, i, G.adjoint());
335 m_S.rightCols(m_n - j - 1).applyOnTheLeft(i - 1, i, G.adjoint());
337 if (!is_negligible(m_S(i, j))) {
340 m_S(i, j) = Scalar(0);
342 if (m_computeQZ) m_Q.applyOnTheRight(i - 1, i, G);
345 if (!numext::is_exactly_zero(m_T.coeff(i, i - 1))) {
347 G.makeGivens(m_T.coeff(i, i), m_T.coeff(i, i - 1), &m_T.coeffRef(i, i));
348 m_T.topRows(i).applyOnTheRight(i - 1, i, G.adjoint());
349 m_T.coeffRef(i, i - 1) = Scalar(0);
351 m_S.applyOnTheRight(i - 1, i, G.adjoint());
353 if (m_computeQZ) m_Z.applyOnTheLeft(i - 1, i, G);
359template <
typename MatrixType>
360template <
typename SparseMatrixType_>
362 m_computeQZ = computeQZ;
364 eigen_assert(m_n == A.cols() &&
"A is not a square matrix");
365 eigen_assert(m_n == B.rows() && m_n == B.cols() &&
"B is not a square matrix or B is not of the same size as A");
366 m_isInitialized =
true;
368 hessenbergTriangularSparse(A, B);
372 reduceHessenbergTriangular();
375template <
typename MatrixType_>
377 Index l = m_n - 1, f;
378 unsigned int local_iter = 0;
381 while (l > 0 && local_iter < m_maxIters) {
382 f = findSmallSubdiagEntry(l);
386 m_S.coeffRef(f, f - 1) = Scalar(0);
391 }
else if (f == l - 1) {
393 reduceDiagonal2x2block(f);
397 Index z = findSmallDiagEntry(f, l);
399 push_down_zero_ST(z, l);
401 do_QZ_step(f, m_n - l - 1);
411template <
typename MatrixType_>
414 S << Scalar(0), Scalar(1), Scalar(1), Scalar(0);
415 Vec2 bprime = S * b.adjoint();
416 JacobiRotation<Scalar> J;
417 J.makeGivens(bprime(0), bprime(1));
419 Z.applyOnTheLeft(0, 1, J);
424template <
typename MatrixType_>
428 const auto a = [p,
this](Index i, Index j) {
return m_S(p + i - 1, p + j - 1); };
429 const auto b = [p,
this](Index i, Index j) {
return m_T(p + i - 1, p + j - 1); };
430 const Index m = m_n - p - q;
433 Scalar W1 = a(m - 1, m - 1) / b(m - 1, m - 1) - a(1, 1) / b(1, 1), W2 = a(m, m) / b(m, m) - a(1, 1) / b(1, 1),
434 W3 = a(m, m - 1) / b(m - 1, m - 1);
436 x = (W1 * W2 - a(m - 1, m) / b(m, m) * W3 + W3 * b(m - 1, m) / b(m, m) * a(1, 1) / b(1, 1)) * b(1, 1) / a(2, 1) +
437 a(1, 2) / b(2, 2) - a(1, 1) / b(1, 1) * b(1, 2) / b(2, 2);
438 y = (a(2, 2) / b(2, 2) - a(1, 1) / b(1, 1)) - a(2, 1) / b(1, 1) * b(1, 2) / b(2, 2) - W1 - W2 +
439 W3 * (b(m - 1, m) / b(m, m));
440 z = a(3, 2) / b(2, 2);
442 const PermutationMatrix<3, 3, int> S3(
Vector3i(2, 0, 1));
443 for (Index k = p; k < p + m - 2; k++) {
448 X.makeHouseholder(ess, tau, beta);
452 m_S.template middleRows<3>(k)
453 .rightCols((std::min)(m_n, m_n - k + 1))
454 .applyHouseholderOnTheLeft(ess, tau, m_ws.data());
455 m_T.template middleRows<3>(k).rightCols(m_n - k).applyHouseholderOnTheLeft(ess, tau, m_ws.data());
456 if (m_computeQZ) m_Q.template middleCols<3>(k).applyHouseholderOnTheRight(ess, std::conj(tau), m_ws.data());
459 Vec3 bprime = (m_T.template block<1, 3>(k + 2, k) * S3).adjoint();
460 bprime.makeHouseholder(ess, tau, beta);
461 m_S.template middleCols<3>(k).topRows((std::min)(k + 4, m_n)).applyOnTheRight(S3);
462 m_S.template middleCols<3>(k)
463 .topRows((std::min)(k + 4, m_n))
464 .applyHouseholderOnTheRight(ess, std::conj(tau), m_ws.data());
465 m_S.template middleCols<3>(k).topRows((std::min)(k + 4, m_n)).applyOnTheRight(S3.transpose());
466 m_T.template middleCols<3>(k).topRows((std::min)(k + 3, m_n)).applyOnTheRight(S3);
467 m_T.template middleCols<3>(k)
468 .topRows((std::min)(k + 3, m_n))
469 .applyHouseholderOnTheRight(ess, std::conj(tau), m_ws.data());
470 m_T.template middleCols<3>(k).topRows((std::min)(k + 3, m_n)).applyOnTheRight(S3.transpose());
472 m_Z.template middleRows<3>(k).applyOnTheLeft(S3.transpose());
473 m_Z.template middleRows<3>(k).applyHouseholderOnTheLeft(ess, tau, m_ws.data());
474 m_Z.template middleRows<3>(k).applyOnTheLeft(S3);
476 Mat2 Zk2 = computeZk2(m_T.template block<1, 2>(k + 1, k));
477 m_S.template middleCols<2>(k).topRows((std::min)(k + 4, m_n)).applyOnTheRight(Zk2);
478 m_T.template middleCols<2>(k).topRows((std::min)(k + 3, m_n)).applyOnTheRight(Zk2);
480 if (m_computeQZ) m_Z.template middleRows<2>(k).applyOnTheLeft(Zk2.adjoint());
490 JacobiRotation<Scalar> J;
492 m_S.template middleRows<2>(p + m - 2).applyOnTheLeft(0, 1, J.adjoint());
493 m_T.template middleRows<2>(p + m - 2).applyOnTheLeft(0, 1, J.adjoint());
495 if (m_computeQZ) m_Q.template middleCols<2>(p + m - 2).applyOnTheRight(0, 1, J);
498 Mat2 Zn1 = computeZk2(m_T.template block<1, 2>(p + m - 1, p + m - 2));
499 m_S.template middleCols<2>(p + m - 2).applyOnTheRight(Zn1);
500 m_T.template middleCols<2>(p + m - 2).applyOnTheRight(Zn1);
502 if (m_computeQZ) m_Z.template middleRows<2>(p + m - 2).applyOnTheLeft(Zn1.adjoint());
506template <
typename MatrixType_>
509 Mat2 Si = m_S.template block<2, 2>(i, i), Ti = m_T.template block<2, 2>(i, i);
510 if (is_negligible(Ti(0, 0)) && !is_negligible(Ti(1, 1))) {
513 m_S.applyOnTheLeft(i, i + 1, G.
adjoint());
514 m_T.applyOnTheLeft(i, i + 1, G.
adjoint());
516 if (m_computeQZ) m_Q.applyOnTheRight(i, i + 1, G);
518 }
else if (!is_negligible(Ti(0, 0)) && is_negligible(Ti(1, 1))) {
520 G.
makeGivens(m_S(i + 1, i + 1), m_S(i + 1, i));
521 m_S.applyOnTheRight(i, i + 1, G.
adjoint());
522 m_T.applyOnTheRight(i, i + 1, G.
adjoint());
523 if (m_computeQZ) m_Z.applyOnTheLeft(i, i + 1, G);
524 }
else if (!is_negligible(Ti(0, 0)) && !is_negligible((Ti(1, 1)))) {
525 Scalar mu = Si(0, 0) / Ti(0, 0);
526 Scalar a12_bar = Si(0, 1) - mu * Ti(0, 1);
527 Scalar a22_bar = Si(1, 1) - mu * Ti(1, 1);
528 Scalar p = Scalar(0.5) * (a22_bar / Ti(1, 1) - Ti(0, 1) * Si(1, 0) / (Ti(0, 0) * Ti(1, 1)));
529 RealScalar sgn_p = p.real() >= RealScalar(0) ? RealScalar(1) : RealScalar(-1);
530 Scalar q = Si(1, 0) * a12_bar / (Ti(0, 0) * Ti(1, 1));
531 Scalar r = p * p + q;
532 Scalar lambda = mu + p + sgn_p * numext::sqrt(r);
533 Mat2 E = Si - lambda * Ti;
535 E.rowwise().norm().maxCoeff(&l);
536 JacobiRotation<Scalar> G;
537 G.makeGivens(E(l, 1), E(l, 0));
538 m_S.applyOnTheRight(i, i + 1, G.adjoint());
539 m_T.applyOnTheRight(i, i + 1, G.adjoint());
541 if (m_computeQZ) m_Z.applyOnTheLeft(i, i + 1, G);
543 Mat2 tildeSi = m_S.template block<2, 2>(i, i), tildeTi = m_T.template block<2, 2>(i, i);
544 Mat2 C = tildeSi.norm() < (lambda * tildeTi).norm() ? tildeSi : lambda * tildeTi;
545 G.makeGivens(C(0, 0), C(1, 0));
546 m_S.applyOnTheLeft(i, i + 1, G.adjoint());
547 m_T.applyOnTheLeft(i, i + 1, G.adjoint());
549 if (m_computeQZ) m_Q.applyOnTheRight(i, i + 1, G);
552 if (!is_negligible(m_S(i + 1, i), m_normOfS * NumTraits<RealScalar>::epsilon())) {
555 m_S(i + 1, i) = Scalar(0);
560template <
typename MatrixType_>
564 JacobiRotation<Scalar> J;
565 for (Index j = k + 1; j <= l; j++) {
567 J.makeGivens(m_T(j - 1, j), m_T(j, j), &m_T.coeffRef(j - 1, j));
568 if (m_n - j - 1 > 0) {
569 m_T.rightCols(m_n - j - 1).applyOnTheLeft(j - 1, j, J.adjoint());
571 m_T.coeffRef(j, j) = Scalar(0);
573 m_S.applyOnTheLeft(j - 1, j, J.adjoint());
575 if (m_computeQZ) m_Q.applyOnTheRight(j - 1, j, J);
579 J.makeGivens(std::conj(m_S(j, j - 1)), std::conj(m_S(j, j - 2)));
580 m_S.applyOnTheRight(j - 1, j - 2, J);
581 m_S(j, j - 2) = Scalar(0);
582 m_T.applyOnTheRight(j - 1, j - 2, J);
583 if (m_computeQZ) m_Z.applyOnTheLeft(j - 1, j - 2, J.adjoint());
589 J.makeGivens(std::conj(m_S(l, l)), std::conj(m_S(l, l - 1)));
590 m_S.topRows(l + 1).applyOnTheRight(l, l - 1, J);
592 if (!is_negligible(m_S(l, l - 1), m_normOfS * NumTraits<Scalar>::epsilon())) {
595 m_S(l, l - 1) = Scalar(0);
597 m_T.topRows(l + 1).applyOnTheRight(l, l - 1, J);
599 if (m_computeQZ) m_Z.applyOnTheLeft(l, l - 1, J.adjoint());
602 if (!is_negligible(m_T(l, l)) || !is_negligible(m_S(l, l - 1))) {
605 m_T(l, l) = Scalar(0);
606 m_S(l, l - 1) = Scalar(0);
611template <
typename MatrixType_>
613 const Index size = m_S.cols();
614 m_normOfS = RealScalar(0);
615 m_normOfT = RealScalar(0);
616 for (Index j = 0; j < size; ++j) {
617 m_normOfS += m_S.col(j).segment(0, (std::min)(size, j + 2)).cwiseAbs().sum();
618 m_normOfT += m_T.row(j).segment(j, size - j).cwiseAbs().sum();
624template <
typename MatrixType_>
628 RealScalar s = numext::abs(m_S.coeff(res - 1, res - 1)) + numext::abs(m_S.coeff(res, res));
629 if (s == Scalar(0)) s = m_normOfS;
630 if (numext::abs(m_S.coeff(res, res - 1)) < NumTraits<RealScalar>::epsilon() * s)
break;
639template <
typename MatrixType_>
643 if (numext::abs(m_T.coeff(res, res)) <= NumTraits<RealScalar>::epsilon() * m_normOfT)
break;
Performs a QZ decomposition of a pair of matrices A, B.
Rotation given by a cosine-sine pair.
Definition Jacobi.h:38
JacobiRotation adjoint() const
Definition Jacobi.h:67
void makeGivens(const Scalar &p, const Scalar &q, Scalar *r=0)
Definition Jacobi.h:152
@ NumericalIssue
Definition Constants.h:442
@ Success
Definition Constants.h:440
@ NoConvergence
Definition Constants.h:444
Matrix< int, 3, 1 > Vector3i
3×1 vector of type int.
Definition Matrix.h:477
Namespace containing all symbols from the Eigen library.
Definition B01_Experimental.dox:1
const int Dynamic
Definition Constants.h:25