RealQZ.h 23 KB

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