RealQZ.h 23 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654
  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. #ifndef EIGEN_REAL_QZ_H
  10. #define EIGEN_REAL_QZ_H
  11. namespace Eigen {
  12. /** \eigenvalues_module \ingroup Eigenvalues_Module
  13. *
  14. *
  15. * \class RealQZ
  16. *
  17. * \brief Performs a real QZ decomposition of a pair of square matrices
  18. *
  19. * \tparam _MatrixType the type of the matrix of which we are computing the
  20. * real QZ decomposition; this is expected to be an instantiation of the
  21. * Matrix class template.
  22. *
  23. * Given a real square matrices A and B, this class computes the real QZ
  24. * decomposition: \f$ A = Q S Z \f$, \f$ B = Q T Z \f$ where Q and Z are
  25. * real orthogonal matrixes, T is upper-triangular matrix, and S is upper
  26. * quasi-triangular matrix. An orthogonal matrix is a matrix whose
  27. * inverse is equal to its transpose, \f$ U^{-1} = U^T \f$. A quasi-triangular
  28. * matrix is a block-triangular matrix whose diagonal consists of 1-by-1
  29. * blocks and 2-by-2 blocks where further reduction is impossible due to
  30. * complex eigenvalues.
  31. *
  32. * The eigenvalues of the pencil \f$ A - z B \f$ can be obtained from
  33. * 1x1 and 2x2 blocks on the diagonals of S and T.
  34. *
  35. * Call the function compute() to compute the real QZ decomposition of a
  36. * given pair of matrices. Alternatively, you can use the
  37. * RealQZ(const MatrixType& B, const MatrixType& B, bool computeQZ)
  38. * constructor which computes the real QZ decomposition at construction
  39. * time. Once the decomposition is computed, you can use the matrixS(),
  40. * matrixT(), matrixQ() and matrixZ() functions to retrieve the matrices
  41. * S, T, Q and Z in the decomposition. If computeQZ==false, some time
  42. * is saved by not computing matrices Q and Z.
  43. *
  44. * Example: \include RealQZ_compute.cpp
  45. * Output: \include RealQZ_compute.out
  46. *
  47. * \note The implementation is based on the algorithm in "Matrix Computations"
  48. * by Gene H. Golub and Charles F. Van Loan, and a paper "An algorithm for
  49. * generalized eigenvalue problems" by C.B.Moler and G.W.Stewart.
  50. *
  51. * \sa class RealSchur, class ComplexSchur, class EigenSolver, class ComplexEigenSolver
  52. */
  53. template<typename _MatrixType> class RealQZ
  54. {
  55. public:
  56. typedef _MatrixType MatrixType;
  57. enum {
  58. RowsAtCompileTime = MatrixType::RowsAtCompileTime,
  59. ColsAtCompileTime = MatrixType::ColsAtCompileTime,
  60. Options = MatrixType::Options,
  61. MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
  62. MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
  63. };
  64. typedef typename MatrixType::Scalar Scalar;
  65. typedef std::complex<typename NumTraits<Scalar>::Real> ComplexScalar;
  66. typedef Eigen::Index Index; ///< \deprecated since Eigen 3.3
  67. typedef Matrix<ComplexScalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> EigenvalueType;
  68. typedef Matrix<Scalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> ColumnVectorType;
  69. /** \brief Default constructor.
  70. *
  71. * \param [in] size Positive integer, size of the matrix whose QZ decomposition will be computed.
  72. *
  73. * The default constructor is useful in cases in which the user intends to
  74. * perform decompositions via compute(). The \p size parameter is only
  75. * used as a hint. It is not an error to give a wrong \p size, but it may
  76. * impair performance.
  77. *
  78. * \sa compute() for an example.
  79. */
  80. explicit RealQZ(Index size = RowsAtCompileTime==Dynamic ? 1 : RowsAtCompileTime) :
  81. m_S(size, size),
  82. m_T(size, size),
  83. m_Q(size, size),
  84. m_Z(size, size),
  85. m_workspace(size*2),
  86. m_maxIters(400),
  87. m_isInitialized(false)
  88. { }
  89. /** \brief Constructor; computes real QZ decomposition of given matrices
  90. *
  91. * \param[in] A Matrix A.
  92. * \param[in] B Matrix B.
  93. * \param[in] computeQZ If false, A and Z are not computed.
  94. *
  95. * This constructor calls compute() to compute the QZ decomposition.
  96. */
  97. RealQZ(const MatrixType& A, const MatrixType& B, bool computeQZ = true) :
  98. m_S(A.rows(),A.cols()),
  99. m_T(A.rows(),A.cols()),
  100. m_Q(A.rows(),A.cols()),
  101. m_Z(A.rows(),A.cols()),
  102. m_workspace(A.rows()*2),
  103. m_maxIters(400),
  104. m_isInitialized(false) {
  105. compute(A, B, computeQZ);
  106. }
  107. /** \brief Returns matrix Q in the QZ decomposition.
  108. *
  109. * \returns A const reference to the matrix Q.
  110. */
  111. const MatrixType& matrixQ() const {
  112. eigen_assert(m_isInitialized && "RealQZ is not initialized.");
  113. eigen_assert(m_computeQZ && "The matrices Q and Z have not been computed during the QZ decomposition.");
  114. return m_Q;
  115. }
  116. /** \brief Returns matrix Z in the QZ decomposition.
  117. *
  118. * \returns A const reference to the matrix Z.
  119. */
  120. const MatrixType& matrixZ() const {
  121. eigen_assert(m_isInitialized && "RealQZ is not initialized.");
  122. eigen_assert(m_computeQZ && "The matrices Q and Z have not been computed during the QZ decomposition.");
  123. return m_Z;
  124. }
  125. /** \brief Returns matrix S in the QZ decomposition.
  126. *
  127. * \returns A const reference to the matrix S.
  128. */
  129. const MatrixType& matrixS() const {
  130. eigen_assert(m_isInitialized && "RealQZ is not initialized.");
  131. return m_S;
  132. }
  133. /** \brief Returns matrix S in the QZ decomposition.
  134. *
  135. * \returns A const reference to the matrix S.
  136. */
  137. const MatrixType& matrixT() const {
  138. eigen_assert(m_isInitialized && "RealQZ is not initialized.");
  139. return m_T;
  140. }
  141. /** \brief Computes QZ decomposition of given matrix.
  142. *
  143. * \param[in] A Matrix A.
  144. * \param[in] B Matrix B.
  145. * \param[in] computeQZ If false, A and Z are not computed.
  146. * \returns Reference to \c *this
  147. */
  148. RealQZ& compute(const MatrixType& A, const MatrixType& B, bool computeQZ = true);
  149. /** \brief Reports whether previous computation was successful.
  150. *
  151. * \returns \c Success if computation was succesful, \c NoConvergence otherwise.
  152. */
  153. ComputationInfo info() const
  154. {
  155. eigen_assert(m_isInitialized && "RealQZ is not initialized.");
  156. return m_info;
  157. }
  158. /** \brief Returns number of performed QR-like iterations.
  159. */
  160. Index iterations() const
  161. {
  162. eigen_assert(m_isInitialized && "RealQZ is not initialized.");
  163. return m_global_iter;
  164. }
  165. /** Sets the maximal number of iterations allowed to converge to one eigenvalue
  166. * or decouple the problem.
  167. */
  168. RealQZ& setMaxIterations(Index maxIters)
  169. {
  170. m_maxIters = maxIters;
  171. return *this;
  172. }
  173. private:
  174. MatrixType m_S, m_T, m_Q, m_Z;
  175. Matrix<Scalar,Dynamic,1> m_workspace;
  176. ComputationInfo m_info;
  177. Index m_maxIters;
  178. bool m_isInitialized;
  179. bool m_computeQZ;
  180. Scalar m_normOfT, m_normOfS;
  181. Index m_global_iter;
  182. typedef Matrix<Scalar,3,1> Vector3s;
  183. typedef Matrix<Scalar,2,1> Vector2s;
  184. typedef Matrix<Scalar,2,2> Matrix2s;
  185. typedef JacobiRotation<Scalar> JRs;
  186. void hessenbergTriangular();
  187. void computeNorms();
  188. Index findSmallSubdiagEntry(Index iu);
  189. Index findSmallDiagEntry(Index f, Index l);
  190. void splitOffTwoRows(Index i);
  191. void pushDownZero(Index z, Index f, Index l);
  192. void step(Index f, Index l, Index iter);
  193. }; // RealQZ
  194. /** \internal Reduces S and T to upper Hessenberg - triangular form */
  195. template<typename MatrixType>
  196. void RealQZ<MatrixType>::hessenbergTriangular()
  197. {
  198. const Index dim = m_S.cols();
  199. // perform QR decomposition of T, overwrite T with R, save Q
  200. HouseholderQR<MatrixType> qrT(m_T);
  201. m_T = qrT.matrixQR();
  202. m_T.template triangularView<StrictlyLower>().setZero();
  203. m_Q = qrT.householderQ();
  204. // overwrite S with Q* S
  205. m_S.applyOnTheLeft(m_Q.adjoint());
  206. // init Z as Identity
  207. if (m_computeQZ)
  208. m_Z = MatrixType::Identity(dim,dim);
  209. // reduce S to upper Hessenberg with Givens rotations
  210. for (Index j=0; j<=dim-3; j++) {
  211. for (Index i=dim-1; i>=j+2; i--) {
  212. JRs G;
  213. // kill S(i,j)
  214. if(m_S.coeff(i,j) != 0)
  215. {
  216. G.makeGivens(m_S.coeff(i-1,j), m_S.coeff(i,j), &m_S.coeffRef(i-1, j));
  217. m_S.coeffRef(i,j) = Scalar(0.0);
  218. m_S.rightCols(dim-j-1).applyOnTheLeft(i-1,i,G.adjoint());
  219. m_T.rightCols(dim-i+1).applyOnTheLeft(i-1,i,G.adjoint());
  220. // update Q
  221. if (m_computeQZ)
  222. m_Q.applyOnTheRight(i-1,i,G);
  223. }
  224. // kill T(i,i-1)
  225. if(m_T.coeff(i,i-1)!=Scalar(0))
  226. {
  227. G.makeGivens(m_T.coeff(i,i), m_T.coeff(i,i-1), &m_T.coeffRef(i,i));
  228. m_T.coeffRef(i,i-1) = Scalar(0.0);
  229. m_S.applyOnTheRight(i,i-1,G);
  230. m_T.topRows(i).applyOnTheRight(i,i-1,G);
  231. // update Z
  232. if (m_computeQZ)
  233. m_Z.applyOnTheLeft(i,i-1,G.adjoint());
  234. }
  235. }
  236. }
  237. }
  238. /** \internal Computes vector L1 norms of S and T when in Hessenberg-Triangular form already */
  239. template<typename MatrixType>
  240. inline void RealQZ<MatrixType>::computeNorms()
  241. {
  242. const Index size = m_S.cols();
  243. m_normOfS = Scalar(0.0);
  244. m_normOfT = Scalar(0.0);
  245. for (Index j = 0; j < size; ++j)
  246. {
  247. m_normOfS += m_S.col(j).segment(0, (std::min)(size,j+2)).cwiseAbs().sum();
  248. m_normOfT += m_T.row(j).segment(j, size - j).cwiseAbs().sum();
  249. }
  250. }
  251. /** \internal Look for single small sub-diagonal element S(res, res-1) and return res (or 0) */
  252. template<typename MatrixType>
  253. inline Index RealQZ<MatrixType>::findSmallSubdiagEntry(Index iu)
  254. {
  255. using std::abs;
  256. Index res = iu;
  257. while (res > 0)
  258. {
  259. Scalar s = abs(m_S.coeff(res-1,res-1)) + abs(m_S.coeff(res,res));
  260. if (s == Scalar(0.0))
  261. s = m_normOfS;
  262. if (abs(m_S.coeff(res,res-1)) < NumTraits<Scalar>::epsilon() * s)
  263. break;
  264. res--;
  265. }
  266. return res;
  267. }
  268. /** \internal Look for single small diagonal element T(res, res) for res between f and l, and return res (or f-1) */
  269. template<typename MatrixType>
  270. inline Index RealQZ<MatrixType>::findSmallDiagEntry(Index f, Index l)
  271. {
  272. using std::abs;
  273. Index res = l;
  274. while (res >= f) {
  275. if (abs(m_T.coeff(res,res)) <= NumTraits<Scalar>::epsilon() * m_normOfT)
  276. break;
  277. res--;
  278. }
  279. return res;
  280. }
  281. /** \internal decouple 2x2 diagonal block in rows i, i+1 if eigenvalues are real */
  282. template<typename MatrixType>
  283. inline void RealQZ<MatrixType>::splitOffTwoRows(Index i)
  284. {
  285. using std::abs;
  286. using std::sqrt;
  287. const Index dim=m_S.cols();
  288. if (abs(m_S.coeff(i+1,i))==Scalar(0))
  289. return;
  290. Index j = findSmallDiagEntry(i,i+1);
  291. if (j==i-1)
  292. {
  293. // block of (S T^{-1})
  294. Matrix2s STi = m_T.template block<2,2>(i,i).template triangularView<Upper>().
  295. template solve<OnTheRight>(m_S.template block<2,2>(i,i));
  296. Scalar p = Scalar(0.5)*(STi(0,0)-STi(1,1));
  297. Scalar q = p*p + STi(1,0)*STi(0,1);
  298. if (q>=0) {
  299. Scalar z = sqrt(q);
  300. // one QR-like iteration for ABi - lambda I
  301. // is enough - when we know exact eigenvalue in advance,
  302. // convergence is immediate
  303. JRs G;
  304. if (p>=0)
  305. G.makeGivens(p + z, STi(1,0));
  306. else
  307. G.makeGivens(p - z, STi(1,0));
  308. m_S.rightCols(dim-i).applyOnTheLeft(i,i+1,G.adjoint());
  309. m_T.rightCols(dim-i).applyOnTheLeft(i,i+1,G.adjoint());
  310. // update Q
  311. if (m_computeQZ)
  312. m_Q.applyOnTheRight(i,i+1,G);
  313. G.makeGivens(m_T.coeff(i+1,i+1), m_T.coeff(i+1,i));
  314. m_S.topRows(i+2).applyOnTheRight(i+1,i,G);
  315. m_T.topRows(i+2).applyOnTheRight(i+1,i,G);
  316. // update Z
  317. if (m_computeQZ)
  318. m_Z.applyOnTheLeft(i+1,i,G.adjoint());
  319. m_S.coeffRef(i+1,i) = Scalar(0.0);
  320. m_T.coeffRef(i+1,i) = Scalar(0.0);
  321. }
  322. }
  323. else
  324. {
  325. pushDownZero(j,i,i+1);
  326. }
  327. }
  328. /** \internal use zero in T(z,z) to zero S(l,l-1), working in block f..l */
  329. template<typename MatrixType>
  330. inline void RealQZ<MatrixType>::pushDownZero(Index z, Index f, Index l)
  331. {
  332. JRs G;
  333. const Index dim = m_S.cols();
  334. for (Index zz=z; zz<l; zz++)
  335. {
  336. // push 0 down
  337. Index firstColS = zz>f ? (zz-1) : zz;
  338. G.makeGivens(m_T.coeff(zz, zz+1), m_T.coeff(zz+1, zz+1));
  339. m_S.rightCols(dim-firstColS).applyOnTheLeft(zz,zz+1,G.adjoint());
  340. m_T.rightCols(dim-zz).applyOnTheLeft(zz,zz+1,G.adjoint());
  341. m_T.coeffRef(zz+1,zz+1) = Scalar(0.0);
  342. // update Q
  343. if (m_computeQZ)
  344. m_Q.applyOnTheRight(zz,zz+1,G);
  345. // kill S(zz+1, zz-1)
  346. if (zz>f)
  347. {
  348. G.makeGivens(m_S.coeff(zz+1, zz), m_S.coeff(zz+1,zz-1));
  349. m_S.topRows(zz+2).applyOnTheRight(zz, zz-1,G);
  350. m_T.topRows(zz+1).applyOnTheRight(zz, zz-1,G);
  351. m_S.coeffRef(zz+1,zz-1) = Scalar(0.0);
  352. // update Z
  353. if (m_computeQZ)
  354. m_Z.applyOnTheLeft(zz,zz-1,G.adjoint());
  355. }
  356. }
  357. // finally kill S(l,l-1)
  358. G.makeGivens(m_S.coeff(l,l), m_S.coeff(l,l-1));
  359. m_S.applyOnTheRight(l,l-1,G);
  360. m_T.applyOnTheRight(l,l-1,G);
  361. m_S.coeffRef(l,l-1)=Scalar(0.0);
  362. // update Z
  363. if (m_computeQZ)
  364. m_Z.applyOnTheLeft(l,l-1,G.adjoint());
  365. }
  366. /** \internal QR-like iterative step for block f..l */
  367. template<typename MatrixType>
  368. inline void RealQZ<MatrixType>::step(Index f, Index l, Index iter)
  369. {
  370. using std::abs;
  371. const Index dim = m_S.cols();
  372. // x, y, z
  373. Scalar x, y, z;
  374. if (iter==10)
  375. {
  376. // Wilkinson ad hoc shift
  377. const Scalar
  378. a11=m_S.coeff(f+0,f+0), a12=m_S.coeff(f+0,f+1),
  379. a21=m_S.coeff(f+1,f+0), a22=m_S.coeff(f+1,f+1), a32=m_S.coeff(f+2,f+1),
  380. b12=m_T.coeff(f+0,f+1),
  381. b11i=Scalar(1.0)/m_T.coeff(f+0,f+0),
  382. b22i=Scalar(1.0)/m_T.coeff(f+1,f+1),
  383. a87=m_S.coeff(l-1,l-2),
  384. a98=m_S.coeff(l-0,l-1),
  385. b77i=Scalar(1.0)/m_T.coeff(l-2,l-2),
  386. b88i=Scalar(1.0)/m_T.coeff(l-1,l-1);
  387. Scalar ss = abs(a87*b77i) + abs(a98*b88i),
  388. lpl = Scalar(1.5)*ss,
  389. ll = ss*ss;
  390. x = ll + a11*a11*b11i*b11i - lpl*a11*b11i + a12*a21*b11i*b22i
  391. - a11*a21*b12*b11i*b11i*b22i;
  392. y = a11*a21*b11i*b11i - lpl*a21*b11i + a21*a22*b11i*b22i
  393. - a21*a21*b12*b11i*b11i*b22i;
  394. z = a21*a32*b11i*b22i;
  395. }
  396. else if (iter==16)
  397. {
  398. // another exceptional shift
  399. 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) /
  400. (m_T.coeff(l-1,l-1)*m_T.coeff(l,l));
  401. y = m_S.coeff(f+1,f)/m_T.coeff(f,f);
  402. z = 0;
  403. }
  404. else if (iter>23 && !(iter%8))
  405. {
  406. // extremely exceptional shift
  407. x = internal::random<Scalar>(-1.0,1.0);
  408. y = internal::random<Scalar>(-1.0,1.0);
  409. z = internal::random<Scalar>(-1.0,1.0);
  410. }
  411. else
  412. {
  413. // Compute the shifts: (x,y,z,0...) = (AB^-1 - l1 I) (AB^-1 - l2 I) e1
  414. // where l1 and l2 are the eigenvalues of the 2x2 matrix C = U V^-1 where
  415. // U and V are 2x2 bottom right sub matrices of A and B. Thus:
  416. // = AB^-1AB^-1 + l1 l2 I - (l1+l2)(AB^-1)
  417. // = AB^-1AB^-1 + det(M) - tr(M)(AB^-1)
  418. // Since we are only interested in having x, y, z with a correct ratio, we have:
  419. const Scalar
  420. a11 = m_S.coeff(f,f), a12 = m_S.coeff(f,f+1),
  421. a21 = m_S.coeff(f+1,f), a22 = m_S.coeff(f+1,f+1),
  422. a32 = m_S.coeff(f+2,f+1),
  423. a88 = m_S.coeff(l-1,l-1), a89 = m_S.coeff(l-1,l),
  424. a98 = m_S.coeff(l,l-1), a99 = m_S.coeff(l,l),
  425. b11 = m_T.coeff(f,f), b12 = m_T.coeff(f,f+1),
  426. b22 = m_T.coeff(f+1,f+1),
  427. b88 = m_T.coeff(l-1,l-1), b89 = m_T.coeff(l-1,l),
  428. b99 = m_T.coeff(l,l);
  429. x = ( (a88/b88 - a11/b11)*(a99/b99 - a11/b11) - (a89/b99)*(a98/b88) + (a98/b88)*(b89/b99)*(a11/b11) ) * (b11/a21)
  430. + a12/b22 - (a11/b11)*(b12/b22);
  431. y = (a22/b22-a11/b11) - (a21/b11)*(b12/b22) - (a88/b88-a11/b11) - (a99/b99-a11/b11) + (a98/b88)*(b89/b99);
  432. z = a32/b22;
  433. }
  434. JRs G;
  435. for (Index k=f; k<=l-2; k++)
  436. {
  437. // variables for Householder reflections
  438. Vector2s essential2;
  439. Scalar tau, beta;
  440. Vector3s hr(x,y,z);
  441. // Q_k to annihilate S(k+1,k-1) and S(k+2,k-1)
  442. hr.makeHouseholderInPlace(tau, beta);
  443. essential2 = hr.template bottomRows<2>();
  444. Index fc=(std::max)(k-1,Index(0)); // first col to update
  445. m_S.template middleRows<3>(k).rightCols(dim-fc).applyHouseholderOnTheLeft(essential2, tau, m_workspace.data());
  446. m_T.template middleRows<3>(k).rightCols(dim-fc).applyHouseholderOnTheLeft(essential2, tau, m_workspace.data());
  447. if (m_computeQZ)
  448. m_Q.template middleCols<3>(k).applyHouseholderOnTheRight(essential2, tau, m_workspace.data());
  449. if (k>f)
  450. m_S.coeffRef(k+2,k-1) = m_S.coeffRef(k+1,k-1) = Scalar(0.0);
  451. // Z_{k1} to annihilate T(k+2,k+1) and T(k+2,k)
  452. hr << m_T.coeff(k+2,k+2),m_T.coeff(k+2,k),m_T.coeff(k+2,k+1);
  453. hr.makeHouseholderInPlace(tau, beta);
  454. essential2 = hr.template bottomRows<2>();
  455. {
  456. Index lr = (std::min)(k+4,dim); // last row to update
  457. Map<Matrix<Scalar,Dynamic,1> > tmp(m_workspace.data(),lr);
  458. // S
  459. tmp = m_S.template middleCols<2>(k).topRows(lr) * essential2;
  460. tmp += m_S.col(k+2).head(lr);
  461. m_S.col(k+2).head(lr) -= tau*tmp;
  462. m_S.template middleCols<2>(k).topRows(lr) -= (tau*tmp) * essential2.adjoint();
  463. // T
  464. tmp = m_T.template middleCols<2>(k).topRows(lr) * essential2;
  465. tmp += m_T.col(k+2).head(lr);
  466. m_T.col(k+2).head(lr) -= tau*tmp;
  467. m_T.template middleCols<2>(k).topRows(lr) -= (tau*tmp) * essential2.adjoint();
  468. }
  469. if (m_computeQZ)
  470. {
  471. // Z
  472. Map<Matrix<Scalar,1,Dynamic> > tmp(m_workspace.data(),dim);
  473. tmp = essential2.adjoint()*(m_Z.template middleRows<2>(k));
  474. tmp += m_Z.row(k+2);
  475. m_Z.row(k+2) -= tau*tmp;
  476. m_Z.template middleRows<2>(k) -= essential2 * (tau*tmp);
  477. }
  478. m_T.coeffRef(k+2,k) = m_T.coeffRef(k+2,k+1) = Scalar(0.0);
  479. // Z_{k2} to annihilate T(k+1,k)
  480. G.makeGivens(m_T.coeff(k+1,k+1), m_T.coeff(k+1,k));
  481. m_S.applyOnTheRight(k+1,k,G);
  482. m_T.applyOnTheRight(k+1,k,G);
  483. // update Z
  484. if (m_computeQZ)
  485. m_Z.applyOnTheLeft(k+1,k,G.adjoint());
  486. m_T.coeffRef(k+1,k) = Scalar(0.0);
  487. // update x,y,z
  488. x = m_S.coeff(k+1,k);
  489. y = m_S.coeff(k+2,k);
  490. if (k < l-2)
  491. z = m_S.coeff(k+3,k);
  492. } // loop over k
  493. // Q_{n-1} to annihilate y = S(l,l-2)
  494. G.makeGivens(x,y);
  495. m_S.applyOnTheLeft(l-1,l,G.adjoint());
  496. m_T.applyOnTheLeft(l-1,l,G.adjoint());
  497. if (m_computeQZ)
  498. m_Q.applyOnTheRight(l-1,l,G);
  499. m_S.coeffRef(l,l-2) = Scalar(0.0);
  500. // Z_{n-1} to annihilate T(l,l-1)
  501. G.makeGivens(m_T.coeff(l,l),m_T.coeff(l,l-1));
  502. m_S.applyOnTheRight(l,l-1,G);
  503. m_T.applyOnTheRight(l,l-1,G);
  504. if (m_computeQZ)
  505. m_Z.applyOnTheLeft(l,l-1,G.adjoint());
  506. m_T.coeffRef(l,l-1) = Scalar(0.0);
  507. }
  508. template<typename MatrixType>
  509. RealQZ<MatrixType>& RealQZ<MatrixType>::compute(const MatrixType& A_in, const MatrixType& B_in, bool computeQZ)
  510. {
  511. const Index dim = A_in.cols();
  512. eigen_assert (A_in.rows()==dim && A_in.cols()==dim
  513. && B_in.rows()==dim && B_in.cols()==dim
  514. && "Need square matrices of the same dimension");
  515. m_isInitialized = true;
  516. m_computeQZ = computeQZ;
  517. m_S = A_in; m_T = B_in;
  518. m_workspace.resize(dim*2);
  519. m_global_iter = 0;
  520. // entrance point: hessenberg triangular decomposition
  521. hessenbergTriangular();
  522. // compute L1 vector norms of T, S into m_normOfS, m_normOfT
  523. computeNorms();
  524. Index l = dim-1,
  525. f,
  526. local_iter = 0;
  527. while (l>0 && local_iter<m_maxIters)
  528. {
  529. f = findSmallSubdiagEntry(l);
  530. // now rows and columns f..l (including) decouple from the rest of the problem
  531. if (f>0) m_S.coeffRef(f,f-1) = Scalar(0.0);
  532. if (f == l) // One root found
  533. {
  534. l--;
  535. local_iter = 0;
  536. }
  537. else if (f == l-1) // Two roots found
  538. {
  539. splitOffTwoRows(f);
  540. l -= 2;
  541. local_iter = 0;
  542. }
  543. else // No convergence yet
  544. {
  545. // if there's zero on diagonal of T, we can isolate an eigenvalue with Givens rotations
  546. Index z = findSmallDiagEntry(f,l);
  547. if (z>=f)
  548. {
  549. // zero found
  550. pushDownZero(z,f,l);
  551. }
  552. else
  553. {
  554. // We are sure now that S.block(f,f, l-f+1,l-f+1) is underuced upper-Hessenberg
  555. // and T.block(f,f, l-f+1,l-f+1) is invertible uper-triangular, which allows to
  556. // apply a QR-like iteration to rows and columns f..l.
  557. step(f,l, local_iter);
  558. local_iter++;
  559. m_global_iter++;
  560. }
  561. }
  562. }
  563. // check if we converged before reaching iterations limit
  564. m_info = (local_iter<m_maxIters) ? Success : NoConvergence;
  565. // For each non triangular 2x2 diagonal block of S,
  566. // reduce the respective 2x2 diagonal block of T to positive diagonal form using 2x2 SVD.
  567. // This step is not mandatory for QZ, but it does help further extraction of eigenvalues/eigenvectors,
  568. // and is in par with Lapack/Matlab QZ.
  569. if(m_info==Success)
  570. {
  571. for(Index i=0; i<dim-1; ++i)
  572. {
  573. if(m_S.coeff(i+1, i) != Scalar(0))
  574. {
  575. JacobiRotation<Scalar> j_left, j_right;
  576. internal::real_2x2_jacobi_svd(m_T, i, i+1, &j_left, &j_right);
  577. // Apply resulting Jacobi rotations
  578. m_S.applyOnTheLeft(i,i+1,j_left);
  579. m_S.applyOnTheRight(i,i+1,j_right);
  580. m_T.applyOnTheLeft(i,i+1,j_left);
  581. m_T.applyOnTheRight(i,i+1,j_right);
  582. m_T(i+1,i) = m_T(i,i+1) = Scalar(0);
  583. if(m_computeQZ) {
  584. m_Q.applyOnTheRight(i,i+1,j_left.transpose());
  585. m_Z.applyOnTheLeft(i,i+1,j_right.transpose());
  586. }
  587. i++;
  588. }
  589. }
  590. }
  591. return *this;
  592. } // end compute
  593. } // end namespace Eigen
  594. #endif //EIGEN_REAL_QZ