RealSchur.h 20 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553
  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. #ifndef EIGEN_REAL_SCHUR_H
  11. #define EIGEN_REAL_SCHUR_H
  12. #include "./HessenbergDecomposition.h"
  13. namespace Eigen {
  14. /** \eigenvalues_module \ingroup Eigenvalues_Module
  15. *
  16. *
  17. * \class RealSchur
  18. *
  19. * \brief Performs a real Schur decomposition of a square matrix
  20. *
  21. * \tparam _MatrixType the type of the matrix of which we are computing the
  22. * real Schur decomposition; this is expected to be an instantiation of the
  23. * Matrix class template.
  24. *
  25. * Given a real square matrix A, this class computes the real Schur
  26. * decomposition: \f$ A = U T U^T \f$ where U is a real orthogonal matrix and
  27. * T is a real quasi-triangular matrix. An orthogonal matrix is a matrix whose
  28. * inverse is equal to its transpose, \f$ U^{-1} = U^T \f$. A quasi-triangular
  29. * matrix is a block-triangular matrix whose diagonal consists of 1-by-1
  30. * blocks and 2-by-2 blocks with complex eigenvalues. The eigenvalues of the
  31. * blocks on the diagonal of T are the same as the eigenvalues of the matrix
  32. * A, and thus the real Schur decomposition is used in EigenSolver to compute
  33. * the eigendecomposition of a matrix.
  34. *
  35. * Call the function compute() to compute the real Schur decomposition of a
  36. * given matrix. Alternatively, you can use the RealSchur(const MatrixType&, bool)
  37. * constructor which computes the real Schur decomposition at construction
  38. * time. Once the decomposition is computed, you can use the matrixU() and
  39. * matrixT() functions to retrieve the matrices U and T in the decomposition.
  40. *
  41. * The documentation of RealSchur(const MatrixType&, bool) contains an example
  42. * of the typical use of this class.
  43. *
  44. * \note The implementation is adapted from
  45. * <a href="http://math.nist.gov/javanumerics/jama/">JAMA</a> (public domain).
  46. * Their code is based on EISPACK.
  47. *
  48. * \sa class ComplexSchur, class EigenSolver, class ComplexEigenSolver
  49. */
  50. template<typename _MatrixType> class RealSchur
  51. {
  52. public:
  53. typedef _MatrixType MatrixType;
  54. enum {
  55. RowsAtCompileTime = MatrixType::RowsAtCompileTime,
  56. ColsAtCompileTime = MatrixType::ColsAtCompileTime,
  57. Options = MatrixType::Options,
  58. MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
  59. MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
  60. };
  61. typedef typename MatrixType::Scalar Scalar;
  62. typedef std::complex<typename NumTraits<Scalar>::Real> ComplexScalar;
  63. typedef Eigen::Index Index; ///< \deprecated since Eigen 3.3
  64. typedef Matrix<ComplexScalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> EigenvalueType;
  65. typedef Matrix<Scalar, ColsAtCompileTime, 1, Options & ~RowMajor, MaxColsAtCompileTime, 1> ColumnVectorType;
  66. /** \brief Default constructor.
  67. *
  68. * \param [in] size Positive integer, size of the matrix whose Schur decomposition will be computed.
  69. *
  70. * The default constructor is useful in cases in which the user intends to
  71. * perform decompositions via compute(). The \p size parameter is only
  72. * used as a hint. It is not an error to give a wrong \p size, but it may
  73. * impair performance.
  74. *
  75. * \sa compute() for an example.
  76. */
  77. explicit RealSchur(Index size = RowsAtCompileTime==Dynamic ? 1 : RowsAtCompileTime)
  78. : m_matT(size, size),
  79. m_matU(size, size),
  80. m_workspaceVector(size),
  81. m_hess(size),
  82. m_isInitialized(false),
  83. m_matUisUptodate(false),
  84. m_maxIters(-1)
  85. { }
  86. /** \brief Constructor; computes real Schur decomposition of given matrix.
  87. *
  88. * \param[in] matrix Square matrix whose Schur decomposition is to be computed.
  89. * \param[in] computeU If true, both T and U are computed; if false, only T is computed.
  90. *
  91. * This constructor calls compute() to compute the Schur decomposition.
  92. *
  93. * Example: \include RealSchur_RealSchur_MatrixType.cpp
  94. * Output: \verbinclude RealSchur_RealSchur_MatrixType.out
  95. */
  96. template<typename InputType>
  97. explicit RealSchur(const EigenBase<InputType>& matrix, bool computeU = true)
  98. : m_matT(matrix.rows(),matrix.cols()),
  99. m_matU(matrix.rows(),matrix.cols()),
  100. m_workspaceVector(matrix.rows()),
  101. m_hess(matrix.rows()),
  102. m_isInitialized(false),
  103. m_matUisUptodate(false),
  104. m_maxIters(-1)
  105. {
  106. compute(matrix.derived(), computeU);
  107. }
  108. /** \brief Returns the orthogonal matrix in the Schur decomposition.
  109. *
  110. * \returns A const reference to the matrix U.
  111. *
  112. * \pre Either the constructor RealSchur(const MatrixType&, bool) or the
  113. * member function compute(const MatrixType&, bool) has been called before
  114. * to compute the Schur decomposition of a matrix, and \p computeU was set
  115. * to true (the default value).
  116. *
  117. * \sa RealSchur(const MatrixType&, bool) for an example
  118. */
  119. const MatrixType& matrixU() const
  120. {
  121. eigen_assert(m_isInitialized && "RealSchur is not initialized.");
  122. eigen_assert(m_matUisUptodate && "The matrix U has not been computed during the RealSchur decomposition.");
  123. return m_matU;
  124. }
  125. /** \brief Returns the quasi-triangular matrix in the Schur decomposition.
  126. *
  127. * \returns A const reference to the matrix T.
  128. *
  129. * \pre Either the constructor RealSchur(const MatrixType&, bool) or the
  130. * member function compute(const MatrixType&, bool) has been called before
  131. * to compute the Schur decomposition of a matrix.
  132. *
  133. * \sa RealSchur(const MatrixType&, bool) for an example
  134. */
  135. const MatrixType& matrixT() const
  136. {
  137. eigen_assert(m_isInitialized && "RealSchur is not initialized.");
  138. return m_matT;
  139. }
  140. /** \brief Computes Schur decomposition of given matrix.
  141. *
  142. * \param[in] matrix Square matrix whose Schur decomposition is to be computed.
  143. * \param[in] computeU If true, both T and U are computed; if false, only T is computed.
  144. * \returns Reference to \c *this
  145. *
  146. * The Schur decomposition is computed by first reducing the matrix to
  147. * Hessenberg form using the class HessenbergDecomposition. The Hessenberg
  148. * matrix is then reduced to triangular form by performing Francis QR
  149. * iterations with implicit double shift. The cost of computing the Schur
  150. * decomposition depends on the number of iterations; as a rough guide, it
  151. * may be taken to be \f$25n^3\f$ flops if \a computeU is true and
  152. * \f$10n^3\f$ flops if \a computeU is false.
  153. *
  154. * Example: \include RealSchur_compute.cpp
  155. * Output: \verbinclude RealSchur_compute.out
  156. *
  157. * \sa compute(const MatrixType&, bool, Index)
  158. */
  159. template<typename InputType>
  160. RealSchur& compute(const EigenBase<InputType>& matrix, bool computeU = true);
  161. /** \brief Computes Schur decomposition of a Hessenberg matrix H = Z T Z^T
  162. * \param[in] matrixH Matrix in Hessenberg form H
  163. * \param[in] matrixQ orthogonal matrix Q that transform a matrix A to H : A = Q H Q^T
  164. * \param computeU Computes the matriX U of the Schur vectors
  165. * \return Reference to \c *this
  166. *
  167. * This routine assumes that the matrix is already reduced in Hessenberg form matrixH
  168. * using either the class HessenbergDecomposition or another mean.
  169. * It computes the upper quasi-triangular matrix T of the Schur decomposition of H
  170. * When computeU is true, this routine computes the matrix U such that
  171. * A = U T U^T = (QZ) T (QZ)^T = Q H Q^T where A is the initial matrix
  172. *
  173. * NOTE Q is referenced if computeU is true; so, if the initial orthogonal matrix
  174. * is not available, the user should give an identity matrix (Q.setIdentity())
  175. *
  176. * \sa compute(const MatrixType&, bool)
  177. */
  178. template<typename HessMatrixType, typename OrthMatrixType>
  179. RealSchur& computeFromHessenberg(const HessMatrixType& matrixH, const OrthMatrixType& matrixQ, bool computeU);
  180. /** \brief Reports whether previous computation was successful.
  181. *
  182. * \returns \c Success if computation was succesful, \c NoConvergence otherwise.
  183. */
  184. ComputationInfo info() const
  185. {
  186. eigen_assert(m_isInitialized && "RealSchur is not initialized.");
  187. return m_info;
  188. }
  189. /** \brief Sets the maximum number of iterations allowed.
  190. *
  191. * If not specified by the user, the maximum number of iterations is m_maxIterationsPerRow times the size
  192. * of the matrix.
  193. */
  194. RealSchur& setMaxIterations(Index maxIters)
  195. {
  196. m_maxIters = maxIters;
  197. return *this;
  198. }
  199. /** \brief Returns the maximum number of iterations. */
  200. Index getMaxIterations()
  201. {
  202. return m_maxIters;
  203. }
  204. /** \brief Maximum number of iterations per row.
  205. *
  206. * If not otherwise specified, the maximum number of iterations is this number times the size of the
  207. * matrix. It is currently set to 40.
  208. */
  209. static const int m_maxIterationsPerRow = 40;
  210. private:
  211. MatrixType m_matT;
  212. MatrixType m_matU;
  213. ColumnVectorType m_workspaceVector;
  214. HessenbergDecomposition<MatrixType> m_hess;
  215. ComputationInfo m_info;
  216. bool m_isInitialized;
  217. bool m_matUisUptodate;
  218. Index m_maxIters;
  219. typedef Matrix<Scalar,3,1> Vector3s;
  220. Scalar computeNormOfT();
  221. Index findSmallSubdiagEntry(Index iu, const Scalar& considerAsZero);
  222. void splitOffTwoRows(Index iu, bool computeU, const Scalar& exshift);
  223. void computeShift(Index iu, Index iter, Scalar& exshift, Vector3s& shiftInfo);
  224. void initFrancisQRStep(Index il, Index iu, const Vector3s& shiftInfo, Index& im, Vector3s& firstHouseholderVector);
  225. void performFrancisQRStep(Index il, Index im, Index iu, bool computeU, const Vector3s& firstHouseholderVector, Scalar* workspace);
  226. };
  227. template<typename MatrixType>
  228. template<typename InputType>
  229. RealSchur<MatrixType>& RealSchur<MatrixType>::compute(const EigenBase<InputType>& matrix, bool computeU)
  230. {
  231. const Scalar considerAsZero = (std::numeric_limits<Scalar>::min)();
  232. eigen_assert(matrix.cols() == matrix.rows());
  233. Index maxIters = m_maxIters;
  234. if (maxIters == -1)
  235. maxIters = m_maxIterationsPerRow * matrix.rows();
  236. Scalar scale = matrix.derived().cwiseAbs().maxCoeff();
  237. if(scale<considerAsZero)
  238. {
  239. m_matT.setZero(matrix.rows(),matrix.cols());
  240. if(computeU)
  241. m_matU.setIdentity(matrix.rows(),matrix.cols());
  242. m_info = Success;
  243. m_isInitialized = true;
  244. m_matUisUptodate = computeU;
  245. return *this;
  246. }
  247. // Step 1. Reduce to Hessenberg form
  248. m_hess.compute(matrix.derived()/scale);
  249. // Step 2. Reduce to real Schur form
  250. computeFromHessenberg(m_hess.matrixH(), m_hess.matrixQ(), computeU);
  251. m_matT *= scale;
  252. return *this;
  253. }
  254. template<typename MatrixType>
  255. template<typename HessMatrixType, typename OrthMatrixType>
  256. RealSchur<MatrixType>& RealSchur<MatrixType>::computeFromHessenberg(const HessMatrixType& matrixH, const OrthMatrixType& matrixQ, bool computeU)
  257. {
  258. using std::abs;
  259. m_matT = matrixH;
  260. if(computeU)
  261. m_matU = matrixQ;
  262. Index maxIters = m_maxIters;
  263. if (maxIters == -1)
  264. maxIters = m_maxIterationsPerRow * matrixH.rows();
  265. m_workspaceVector.resize(m_matT.cols());
  266. Scalar* workspace = &m_workspaceVector.coeffRef(0);
  267. // The matrix m_matT is divided in three parts.
  268. // Rows 0,...,il-1 are decoupled from the rest because m_matT(il,il-1) is zero.
  269. // Rows il,...,iu is the part we are working on (the active window).
  270. // Rows iu+1,...,end are already brought in triangular form.
  271. Index iu = m_matT.cols() - 1;
  272. Index iter = 0; // iteration count for current eigenvalue
  273. Index totalIter = 0; // iteration count for whole matrix
  274. Scalar exshift(0); // sum of exceptional shifts
  275. Scalar norm = computeNormOfT();
  276. // sub-diagonal entries smaller than considerAsZero will be treated as zero.
  277. // We use eps^2 to enable more precision in small eigenvalues.
  278. Scalar considerAsZero = numext::maxi<Scalar>( norm * numext::abs2(NumTraits<Scalar>::epsilon()),
  279. (std::numeric_limits<Scalar>::min)() );
  280. if(norm!=Scalar(0))
  281. {
  282. while (iu >= 0)
  283. {
  284. Index il = findSmallSubdiagEntry(iu,considerAsZero);
  285. // Check for convergence
  286. if (il == iu) // One root found
  287. {
  288. m_matT.coeffRef(iu,iu) = m_matT.coeff(iu,iu) + exshift;
  289. if (iu > 0)
  290. m_matT.coeffRef(iu, iu-1) = Scalar(0);
  291. iu--;
  292. iter = 0;
  293. }
  294. else if (il == iu-1) // Two roots found
  295. {
  296. splitOffTwoRows(iu, computeU, exshift);
  297. iu -= 2;
  298. iter = 0;
  299. }
  300. else // No convergence yet
  301. {
  302. // The firstHouseholderVector vector has to be initialized to something to get rid of a silly GCC warning (-O1 -Wall -DNDEBUG )
  303. Vector3s firstHouseholderVector = Vector3s::Zero(), shiftInfo;
  304. computeShift(iu, iter, exshift, shiftInfo);
  305. iter = iter + 1;
  306. totalIter = totalIter + 1;
  307. if (totalIter > maxIters) break;
  308. Index im;
  309. initFrancisQRStep(il, iu, shiftInfo, im, firstHouseholderVector);
  310. performFrancisQRStep(il, im, iu, computeU, firstHouseholderVector, workspace);
  311. }
  312. }
  313. }
  314. if(totalIter <= maxIters)
  315. m_info = Success;
  316. else
  317. m_info = NoConvergence;
  318. m_isInitialized = true;
  319. m_matUisUptodate = computeU;
  320. return *this;
  321. }
  322. /** \internal Computes and returns vector L1 norm of T */
  323. template<typename MatrixType>
  324. inline typename MatrixType::Scalar RealSchur<MatrixType>::computeNormOfT()
  325. {
  326. const Index size = m_matT.cols();
  327. // FIXME to be efficient the following would requires a triangular reduxion code
  328. // Scalar norm = m_matT.upper().cwiseAbs().sum()
  329. // + m_matT.bottomLeftCorner(size-1,size-1).diagonal().cwiseAbs().sum();
  330. Scalar norm(0);
  331. for (Index j = 0; j < size; ++j)
  332. norm += m_matT.col(j).segment(0, (std::min)(size,j+2)).cwiseAbs().sum();
  333. return norm;
  334. }
  335. /** \internal Look for single small sub-diagonal element and returns its index */
  336. template<typename MatrixType>
  337. inline Index RealSchur<MatrixType>::findSmallSubdiagEntry(Index iu, const Scalar& considerAsZero)
  338. {
  339. using std::abs;
  340. Index res = iu;
  341. while (res > 0)
  342. {
  343. Scalar s = abs(m_matT.coeff(res-1,res-1)) + abs(m_matT.coeff(res,res));
  344. s = numext::maxi<Scalar>(s * NumTraits<Scalar>::epsilon(), considerAsZero);
  345. if (abs(m_matT.coeff(res,res-1)) <= s)
  346. break;
  347. res--;
  348. }
  349. return res;
  350. }
  351. /** \internal Update T given that rows iu-1 and iu decouple from the rest. */
  352. template<typename MatrixType>
  353. inline void RealSchur<MatrixType>::splitOffTwoRows(Index iu, bool computeU, const Scalar& exshift)
  354. {
  355. using std::sqrt;
  356. using std::abs;
  357. const Index size = m_matT.cols();
  358. // The eigenvalues of the 2x2 matrix [a b; c d] are
  359. // trace +/- sqrt(discr/4) where discr = tr^2 - 4*det, tr = a + d, det = ad - bc
  360. Scalar p = Scalar(0.5) * (m_matT.coeff(iu-1,iu-1) - m_matT.coeff(iu,iu));
  361. Scalar q = p * p + m_matT.coeff(iu,iu-1) * m_matT.coeff(iu-1,iu); // q = tr^2 / 4 - det = discr/4
  362. m_matT.coeffRef(iu,iu) += exshift;
  363. m_matT.coeffRef(iu-1,iu-1) += exshift;
  364. if (q >= Scalar(0)) // Two real eigenvalues
  365. {
  366. Scalar z = sqrt(abs(q));
  367. JacobiRotation<Scalar> rot;
  368. if (p >= Scalar(0))
  369. rot.makeGivens(p + z, m_matT.coeff(iu, iu-1));
  370. else
  371. rot.makeGivens(p - z, m_matT.coeff(iu, iu-1));
  372. m_matT.rightCols(size-iu+1).applyOnTheLeft(iu-1, iu, rot.adjoint());
  373. m_matT.topRows(iu+1).applyOnTheRight(iu-1, iu, rot);
  374. m_matT.coeffRef(iu, iu-1) = Scalar(0);
  375. if (computeU)
  376. m_matU.applyOnTheRight(iu-1, iu, rot);
  377. }
  378. if (iu > 1)
  379. m_matT.coeffRef(iu-1, iu-2) = Scalar(0);
  380. }
  381. /** \internal Form shift in shiftInfo, and update exshift if an exceptional shift is performed. */
  382. template<typename MatrixType>
  383. inline void RealSchur<MatrixType>::computeShift(Index iu, Index iter, Scalar& exshift, Vector3s& shiftInfo)
  384. {
  385. using std::sqrt;
  386. using std::abs;
  387. shiftInfo.coeffRef(0) = m_matT.coeff(iu,iu);
  388. shiftInfo.coeffRef(1) = m_matT.coeff(iu-1,iu-1);
  389. shiftInfo.coeffRef(2) = m_matT.coeff(iu,iu-1) * m_matT.coeff(iu-1,iu);
  390. // Wilkinson's original ad hoc shift
  391. if (iter == 10)
  392. {
  393. exshift += shiftInfo.coeff(0);
  394. for (Index i = 0; i <= iu; ++i)
  395. m_matT.coeffRef(i,i) -= shiftInfo.coeff(0);
  396. Scalar s = abs(m_matT.coeff(iu,iu-1)) + abs(m_matT.coeff(iu-1,iu-2));
  397. shiftInfo.coeffRef(0) = Scalar(0.75) * s;
  398. shiftInfo.coeffRef(1) = Scalar(0.75) * s;
  399. shiftInfo.coeffRef(2) = Scalar(-0.4375) * s * s;
  400. }
  401. // MATLAB's new ad hoc shift
  402. if (iter == 30)
  403. {
  404. Scalar s = (shiftInfo.coeff(1) - shiftInfo.coeff(0)) / Scalar(2.0);
  405. s = s * s + shiftInfo.coeff(2);
  406. if (s > Scalar(0))
  407. {
  408. s = sqrt(s);
  409. if (shiftInfo.coeff(1) < shiftInfo.coeff(0))
  410. s = -s;
  411. s = s + (shiftInfo.coeff(1) - shiftInfo.coeff(0)) / Scalar(2.0);
  412. s = shiftInfo.coeff(0) - shiftInfo.coeff(2) / s;
  413. exshift += s;
  414. for (Index i = 0; i <= iu; ++i)
  415. m_matT.coeffRef(i,i) -= s;
  416. shiftInfo.setConstant(Scalar(0.964));
  417. }
  418. }
  419. }
  420. /** \internal Compute index im at which Francis QR step starts and the first Householder vector. */
  421. template<typename MatrixType>
  422. inline void RealSchur<MatrixType>::initFrancisQRStep(Index il, Index iu, const Vector3s& shiftInfo, Index& im, Vector3s& firstHouseholderVector)
  423. {
  424. using std::abs;
  425. Vector3s& v = firstHouseholderVector; // alias to save typing
  426. for (im = iu-2; im >= il; --im)
  427. {
  428. const Scalar Tmm = m_matT.coeff(im,im);
  429. const Scalar r = shiftInfo.coeff(0) - Tmm;
  430. const Scalar s = shiftInfo.coeff(1) - Tmm;
  431. v.coeffRef(0) = (r * s - shiftInfo.coeff(2)) / m_matT.coeff(im+1,im) + m_matT.coeff(im,im+1);
  432. v.coeffRef(1) = m_matT.coeff(im+1,im+1) - Tmm - r - s;
  433. v.coeffRef(2) = m_matT.coeff(im+2,im+1);
  434. if (im == il) {
  435. break;
  436. }
  437. const Scalar lhs = m_matT.coeff(im,im-1) * (abs(v.coeff(1)) + abs(v.coeff(2)));
  438. const Scalar rhs = v.coeff(0) * (abs(m_matT.coeff(im-1,im-1)) + abs(Tmm) + abs(m_matT.coeff(im+1,im+1)));
  439. if (abs(lhs) < NumTraits<Scalar>::epsilon() * rhs)
  440. break;
  441. }
  442. }
  443. /** \internal Perform a Francis QR step involving rows il:iu and columns im:iu. */
  444. template<typename MatrixType>
  445. inline void RealSchur<MatrixType>::performFrancisQRStep(Index il, Index im, Index iu, bool computeU, const Vector3s& firstHouseholderVector, Scalar* workspace)
  446. {
  447. eigen_assert(im >= il);
  448. eigen_assert(im <= iu-2);
  449. const Index size = m_matT.cols();
  450. for (Index k = im; k <= iu-2; ++k)
  451. {
  452. bool firstIteration = (k == im);
  453. Vector3s v;
  454. if (firstIteration)
  455. v = firstHouseholderVector;
  456. else
  457. v = m_matT.template block<3,1>(k,k-1);
  458. Scalar tau, beta;
  459. Matrix<Scalar, 2, 1> ess;
  460. v.makeHouseholder(ess, tau, beta);
  461. if (beta != Scalar(0)) // if v is not zero
  462. {
  463. if (firstIteration && k > il)
  464. m_matT.coeffRef(k,k-1) = -m_matT.coeff(k,k-1);
  465. else if (!firstIteration)
  466. m_matT.coeffRef(k,k-1) = beta;
  467. // These Householder transformations form the O(n^3) part of the algorithm
  468. m_matT.block(k, k, 3, size-k).applyHouseholderOnTheLeft(ess, tau, workspace);
  469. m_matT.block(0, k, (std::min)(iu,k+3) + 1, 3).applyHouseholderOnTheRight(ess, tau, workspace);
  470. if (computeU)
  471. m_matU.block(0, k, size, 3).applyHouseholderOnTheRight(ess, tau, workspace);
  472. }
  473. }
  474. Matrix<Scalar, 2, 1> v = m_matT.template block<2,1>(iu-1, iu-2);
  475. Scalar tau, beta;
  476. Matrix<Scalar, 1, 1> ess;
  477. v.makeHouseholder(ess, tau, beta);
  478. if (beta != Scalar(0)) // if v is not zero
  479. {
  480. m_matT.coeffRef(iu-1, iu-2) = beta;
  481. m_matT.block(iu-1, iu-1, 2, size-iu+1).applyHouseholderOnTheLeft(ess, tau, workspace);
  482. m_matT.block(0, iu-1, iu+1, 2).applyHouseholderOnTheRight(ess, tau, workspace);
  483. if (computeU)
  484. m_matU.block(0, iu-1, size, 2).applyHouseholderOnTheRight(ess, tau, workspace);
  485. }
  486. // clean up pollution due to round-off errors
  487. for (Index i = im+2; i <= iu; ++i)
  488. {
  489. m_matT.coeffRef(i,i-2) = Scalar(0);
  490. if (i > im+2)
  491. m_matT.coeffRef(i,i-3) = Scalar(0);
  492. }
  493. }
  494. } // end namespace Eigen
  495. #endif // EIGEN_REAL_SCHUR_H