1*bf2c3715SXin Li // This file is part of Eigen, a lightweight C++ template library 2*bf2c3715SXin Li // for linear algebra. 3*bf2c3715SXin Li // 4*bf2c3715SXin Li // Copyright (C) 2009 Claire Maurice 5*bf2c3715SXin Li // Copyright (C) 2009 Gael Guennebaud <[email protected]> 6*bf2c3715SXin Li // Copyright (C) 2010,2012 Jitse Niesen <[email protected]> 7*bf2c3715SXin Li // 8*bf2c3715SXin Li // This Source Code Form is subject to the terms of the Mozilla 9*bf2c3715SXin Li // Public License v. 2.0. If a copy of the MPL was not distributed 10*bf2c3715SXin Li // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. 11*bf2c3715SXin Li 12*bf2c3715SXin Li #ifndef EIGEN_COMPLEX_SCHUR_H 13*bf2c3715SXin Li #define EIGEN_COMPLEX_SCHUR_H 14*bf2c3715SXin Li 15*bf2c3715SXin Li #include "./HessenbergDecomposition.h" 16*bf2c3715SXin Li 17*bf2c3715SXin Li namespace Eigen { 18*bf2c3715SXin Li 19*bf2c3715SXin Li namespace internal { 20*bf2c3715SXin Li template<typename MatrixType, bool IsComplex> struct complex_schur_reduce_to_hessenberg; 21*bf2c3715SXin Li } 22*bf2c3715SXin Li 23*bf2c3715SXin Li /** \eigenvalues_module \ingroup Eigenvalues_Module 24*bf2c3715SXin Li * 25*bf2c3715SXin Li * 26*bf2c3715SXin Li * \class ComplexSchur 27*bf2c3715SXin Li * 28*bf2c3715SXin Li * \brief Performs a complex Schur decomposition of a real or complex square matrix 29*bf2c3715SXin Li * 30*bf2c3715SXin Li * \tparam _MatrixType the type of the matrix of which we are 31*bf2c3715SXin Li * computing the Schur decomposition; this is expected to be an 32*bf2c3715SXin Li * instantiation of the Matrix class template. 33*bf2c3715SXin Li * 34*bf2c3715SXin Li * Given a real or complex square matrix A, this class computes the 35*bf2c3715SXin Li * Schur decomposition: \f$ A = U T U^*\f$ where U is a unitary 36*bf2c3715SXin Li * complex matrix, and T is a complex upper triangular matrix. The 37*bf2c3715SXin Li * diagonal of the matrix T corresponds to the eigenvalues of the 38*bf2c3715SXin Li * matrix A. 39*bf2c3715SXin Li * 40*bf2c3715SXin Li * Call the function compute() to compute the Schur decomposition of 41*bf2c3715SXin Li * a given matrix. Alternatively, you can use the 42*bf2c3715SXin Li * ComplexSchur(const MatrixType&, bool) constructor which computes 43*bf2c3715SXin Li * the Schur decomposition at construction time. Once the 44*bf2c3715SXin Li * decomposition is computed, you can use the matrixU() and matrixT() 45*bf2c3715SXin Li * functions to retrieve the matrices U and V in the decomposition. 46*bf2c3715SXin Li * 47*bf2c3715SXin Li * \note This code is inspired from Jampack 48*bf2c3715SXin Li * 49*bf2c3715SXin Li * \sa class RealSchur, class EigenSolver, class ComplexEigenSolver 50*bf2c3715SXin Li */ 51*bf2c3715SXin Li template<typename _MatrixType> class ComplexSchur 52*bf2c3715SXin Li { 53*bf2c3715SXin Li public: 54*bf2c3715SXin Li typedef _MatrixType MatrixType; 55*bf2c3715SXin Li enum { 56*bf2c3715SXin Li RowsAtCompileTime = MatrixType::RowsAtCompileTime, 57*bf2c3715SXin Li ColsAtCompileTime = MatrixType::ColsAtCompileTime, 58*bf2c3715SXin Li Options = MatrixType::Options, 59*bf2c3715SXin Li MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime, 60*bf2c3715SXin Li MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime 61*bf2c3715SXin Li }; 62*bf2c3715SXin Li 63*bf2c3715SXin Li /** \brief Scalar type for matrices of type \p _MatrixType. */ 64*bf2c3715SXin Li typedef typename MatrixType::Scalar Scalar; 65*bf2c3715SXin Li typedef typename NumTraits<Scalar>::Real RealScalar; 66*bf2c3715SXin Li typedef Eigen::Index Index; ///< \deprecated since Eigen 3.3 67*bf2c3715SXin Li 68*bf2c3715SXin Li /** \brief Complex scalar type for \p _MatrixType. 69*bf2c3715SXin Li * 70*bf2c3715SXin Li * This is \c std::complex<Scalar> if #Scalar is real (e.g., 71*bf2c3715SXin Li * \c float or \c double) and just \c Scalar if #Scalar is 72*bf2c3715SXin Li * complex. 73*bf2c3715SXin Li */ 74*bf2c3715SXin Li typedef std::complex<RealScalar> ComplexScalar; 75*bf2c3715SXin Li 76*bf2c3715SXin Li /** \brief Type for the matrices in the Schur decomposition. 77*bf2c3715SXin Li * 78*bf2c3715SXin Li * This is a square matrix with entries of type #ComplexScalar. 79*bf2c3715SXin Li * The size is the same as the size of \p _MatrixType. 80*bf2c3715SXin Li */ 81*bf2c3715SXin Li typedef Matrix<ComplexScalar, RowsAtCompileTime, ColsAtCompileTime, Options, MaxRowsAtCompileTime, MaxColsAtCompileTime> ComplexMatrixType; 82*bf2c3715SXin Li 83*bf2c3715SXin Li /** \brief Default constructor. 84*bf2c3715SXin Li * 85*bf2c3715SXin Li * \param [in] size Positive integer, size of the matrix whose Schur decomposition will be computed. 86*bf2c3715SXin Li * 87*bf2c3715SXin Li * The default constructor is useful in cases in which the user 88*bf2c3715SXin Li * intends to perform decompositions via compute(). The \p size 89*bf2c3715SXin Li * parameter is only used as a hint. It is not an error to give a 90*bf2c3715SXin Li * wrong \p size, but it may impair performance. 91*bf2c3715SXin Li * 92*bf2c3715SXin Li * \sa compute() for an example. 93*bf2c3715SXin Li */ 94*bf2c3715SXin Li explicit ComplexSchur(Index size = RowsAtCompileTime==Dynamic ? 1 : RowsAtCompileTime) m_matT(size,size)95*bf2c3715SXin Li : m_matT(size,size), 96*bf2c3715SXin Li m_matU(size,size), 97*bf2c3715SXin Li m_hess(size), 98*bf2c3715SXin Li m_isInitialized(false), 99*bf2c3715SXin Li m_matUisUptodate(false), 100*bf2c3715SXin Li m_maxIters(-1) 101*bf2c3715SXin Li {} 102*bf2c3715SXin Li 103*bf2c3715SXin Li /** \brief Constructor; computes Schur decomposition of given matrix. 104*bf2c3715SXin Li * 105*bf2c3715SXin Li * \param[in] matrix Square matrix whose Schur decomposition is to be computed. 106*bf2c3715SXin Li * \param[in] computeU If true, both T and U are computed; if false, only T is computed. 107*bf2c3715SXin Li * 108*bf2c3715SXin Li * This constructor calls compute() to compute the Schur decomposition. 109*bf2c3715SXin Li * 110*bf2c3715SXin Li * \sa matrixT() and matrixU() for examples. 111*bf2c3715SXin Li */ 112*bf2c3715SXin Li template<typename InputType> 113*bf2c3715SXin Li explicit ComplexSchur(const EigenBase<InputType>& matrix, bool computeU = true) 114*bf2c3715SXin Li : m_matT(matrix.rows(),matrix.cols()), 115*bf2c3715SXin Li m_matU(matrix.rows(),matrix.cols()), 116*bf2c3715SXin Li m_hess(matrix.rows()), 117*bf2c3715SXin Li m_isInitialized(false), 118*bf2c3715SXin Li m_matUisUptodate(false), 119*bf2c3715SXin Li m_maxIters(-1) 120*bf2c3715SXin Li { 121*bf2c3715SXin Li compute(matrix.derived(), computeU); 122*bf2c3715SXin Li } 123*bf2c3715SXin Li 124*bf2c3715SXin Li /** \brief Returns the unitary matrix in the Schur decomposition. 125*bf2c3715SXin Li * 126*bf2c3715SXin Li * \returns A const reference to the matrix U. 127*bf2c3715SXin Li * 128*bf2c3715SXin Li * It is assumed that either the constructor 129*bf2c3715SXin Li * ComplexSchur(const MatrixType& matrix, bool computeU) or the 130*bf2c3715SXin Li * member function compute(const MatrixType& matrix, bool computeU) 131*bf2c3715SXin Li * has been called before to compute the Schur decomposition of a 132*bf2c3715SXin Li * matrix, and that \p computeU was set to true (the default 133*bf2c3715SXin Li * value). 134*bf2c3715SXin Li * 135*bf2c3715SXin Li * Example: \include ComplexSchur_matrixU.cpp 136*bf2c3715SXin Li * Output: \verbinclude ComplexSchur_matrixU.out 137*bf2c3715SXin Li */ matrixU()138*bf2c3715SXin Li const ComplexMatrixType& matrixU() const 139*bf2c3715SXin Li { 140*bf2c3715SXin Li eigen_assert(m_isInitialized && "ComplexSchur is not initialized."); 141*bf2c3715SXin Li eigen_assert(m_matUisUptodate && "The matrix U has not been computed during the ComplexSchur decomposition."); 142*bf2c3715SXin Li return m_matU; 143*bf2c3715SXin Li } 144*bf2c3715SXin Li 145*bf2c3715SXin Li /** \brief Returns the triangular matrix in the Schur decomposition. 146*bf2c3715SXin Li * 147*bf2c3715SXin Li * \returns A const reference to the matrix T. 148*bf2c3715SXin Li * 149*bf2c3715SXin Li * It is assumed that either the constructor 150*bf2c3715SXin Li * ComplexSchur(const MatrixType& matrix, bool computeU) or the 151*bf2c3715SXin Li * member function compute(const MatrixType& matrix, bool computeU) 152*bf2c3715SXin Li * has been called before to compute the Schur decomposition of a 153*bf2c3715SXin Li * matrix. 154*bf2c3715SXin Li * 155*bf2c3715SXin Li * Note that this function returns a plain square matrix. If you want to reference 156*bf2c3715SXin Li * only the upper triangular part, use: 157*bf2c3715SXin Li * \code schur.matrixT().triangularView<Upper>() \endcode 158*bf2c3715SXin Li * 159*bf2c3715SXin Li * Example: \include ComplexSchur_matrixT.cpp 160*bf2c3715SXin Li * Output: \verbinclude ComplexSchur_matrixT.out 161*bf2c3715SXin Li */ matrixT()162*bf2c3715SXin Li const ComplexMatrixType& matrixT() const 163*bf2c3715SXin Li { 164*bf2c3715SXin Li eigen_assert(m_isInitialized && "ComplexSchur is not initialized."); 165*bf2c3715SXin Li return m_matT; 166*bf2c3715SXin Li } 167*bf2c3715SXin Li 168*bf2c3715SXin Li /** \brief Computes Schur decomposition of given matrix. 169*bf2c3715SXin Li * 170*bf2c3715SXin Li * \param[in] matrix Square matrix whose Schur decomposition is to be computed. 171*bf2c3715SXin Li * \param[in] computeU If true, both T and U are computed; if false, only T is computed. 172*bf2c3715SXin Li 173*bf2c3715SXin Li * \returns Reference to \c *this 174*bf2c3715SXin Li * 175*bf2c3715SXin Li * The Schur decomposition is computed by first reducing the 176*bf2c3715SXin Li * matrix to Hessenberg form using the class 177*bf2c3715SXin Li * HessenbergDecomposition. The Hessenberg matrix is then reduced 178*bf2c3715SXin Li * to triangular form by performing QR iterations with a single 179*bf2c3715SXin Li * shift. The cost of computing the Schur decomposition depends 180*bf2c3715SXin Li * on the number of iterations; as a rough guide, it may be taken 181*bf2c3715SXin Li * on the number of iterations; as a rough guide, it may be taken 182*bf2c3715SXin Li * to be \f$25n^3\f$ complex flops, or \f$10n^3\f$ complex flops 183*bf2c3715SXin Li * if \a computeU is false. 184*bf2c3715SXin Li * 185*bf2c3715SXin Li * Example: \include ComplexSchur_compute.cpp 186*bf2c3715SXin Li * Output: \verbinclude ComplexSchur_compute.out 187*bf2c3715SXin Li * 188*bf2c3715SXin Li * \sa compute(const MatrixType&, bool, Index) 189*bf2c3715SXin Li */ 190*bf2c3715SXin Li template<typename InputType> 191*bf2c3715SXin Li ComplexSchur& compute(const EigenBase<InputType>& matrix, bool computeU = true); 192*bf2c3715SXin Li 193*bf2c3715SXin Li /** \brief Compute Schur decomposition from a given Hessenberg matrix 194*bf2c3715SXin Li * \param[in] matrixH Matrix in Hessenberg form H 195*bf2c3715SXin Li * \param[in] matrixQ orthogonal matrix Q that transform a matrix A to H : A = Q H Q^T 196*bf2c3715SXin Li * \param computeU Computes the matriX U of the Schur vectors 197*bf2c3715SXin Li * \return Reference to \c *this 198*bf2c3715SXin Li * 199*bf2c3715SXin Li * This routine assumes that the matrix is already reduced in Hessenberg form matrixH 200*bf2c3715SXin Li * using either the class HessenbergDecomposition or another mean. 201*bf2c3715SXin Li * It computes the upper quasi-triangular matrix T of the Schur decomposition of H 202*bf2c3715SXin Li * When computeU is true, this routine computes the matrix U such that 203*bf2c3715SXin Li * A = U T U^T = (QZ) T (QZ)^T = Q H Q^T where A is the initial matrix 204*bf2c3715SXin Li * 205*bf2c3715SXin Li * NOTE Q is referenced if computeU is true; so, if the initial orthogonal matrix 206*bf2c3715SXin Li * is not available, the user should give an identity matrix (Q.setIdentity()) 207*bf2c3715SXin Li * 208*bf2c3715SXin Li * \sa compute(const MatrixType&, bool) 209*bf2c3715SXin Li */ 210*bf2c3715SXin Li template<typename HessMatrixType, typename OrthMatrixType> 211*bf2c3715SXin Li ComplexSchur& computeFromHessenberg(const HessMatrixType& matrixH, const OrthMatrixType& matrixQ, bool computeU=true); 212*bf2c3715SXin Li 213*bf2c3715SXin Li /** \brief Reports whether previous computation was successful. 214*bf2c3715SXin Li * 215*bf2c3715SXin Li * \returns \c Success if computation was successful, \c NoConvergence otherwise. 216*bf2c3715SXin Li */ info()217*bf2c3715SXin Li ComputationInfo info() const 218*bf2c3715SXin Li { 219*bf2c3715SXin Li eigen_assert(m_isInitialized && "ComplexSchur is not initialized."); 220*bf2c3715SXin Li return m_info; 221*bf2c3715SXin Li } 222*bf2c3715SXin Li 223*bf2c3715SXin Li /** \brief Sets the maximum number of iterations allowed. 224*bf2c3715SXin Li * 225*bf2c3715SXin Li * If not specified by the user, the maximum number of iterations is m_maxIterationsPerRow times the size 226*bf2c3715SXin Li * of the matrix. 227*bf2c3715SXin Li */ setMaxIterations(Index maxIters)228*bf2c3715SXin Li ComplexSchur& setMaxIterations(Index maxIters) 229*bf2c3715SXin Li { 230*bf2c3715SXin Li m_maxIters = maxIters; 231*bf2c3715SXin Li return *this; 232*bf2c3715SXin Li } 233*bf2c3715SXin Li 234*bf2c3715SXin Li /** \brief Returns the maximum number of iterations. */ getMaxIterations()235*bf2c3715SXin Li Index getMaxIterations() 236*bf2c3715SXin Li { 237*bf2c3715SXin Li return m_maxIters; 238*bf2c3715SXin Li } 239*bf2c3715SXin Li 240*bf2c3715SXin Li /** \brief Maximum number of iterations per row. 241*bf2c3715SXin Li * 242*bf2c3715SXin Li * If not otherwise specified, the maximum number of iterations is this number times the size of the 243*bf2c3715SXin Li * matrix. It is currently set to 30. 244*bf2c3715SXin Li */ 245*bf2c3715SXin Li static const int m_maxIterationsPerRow = 30; 246*bf2c3715SXin Li 247*bf2c3715SXin Li protected: 248*bf2c3715SXin Li ComplexMatrixType m_matT, m_matU; 249*bf2c3715SXin Li HessenbergDecomposition<MatrixType> m_hess; 250*bf2c3715SXin Li ComputationInfo m_info; 251*bf2c3715SXin Li bool m_isInitialized; 252*bf2c3715SXin Li bool m_matUisUptodate; 253*bf2c3715SXin Li Index m_maxIters; 254*bf2c3715SXin Li 255*bf2c3715SXin Li private: 256*bf2c3715SXin Li bool subdiagonalEntryIsNeglegible(Index i); 257*bf2c3715SXin Li ComplexScalar computeShift(Index iu, Index iter); 258*bf2c3715SXin Li void reduceToTriangularForm(bool computeU); 259*bf2c3715SXin Li friend struct internal::complex_schur_reduce_to_hessenberg<MatrixType, NumTraits<Scalar>::IsComplex>; 260*bf2c3715SXin Li }; 261*bf2c3715SXin Li 262*bf2c3715SXin Li /** If m_matT(i+1,i) is neglegible in floating point arithmetic 263*bf2c3715SXin Li * compared to m_matT(i,i) and m_matT(j,j), then set it to zero and 264*bf2c3715SXin Li * return true, else return false. */ 265*bf2c3715SXin Li template<typename MatrixType> 266*bf2c3715SXin Li inline bool ComplexSchur<MatrixType>::subdiagonalEntryIsNeglegible(Index i) 267*bf2c3715SXin Li { 268*bf2c3715SXin Li RealScalar d = numext::norm1(m_matT.coeff(i,i)) + numext::norm1(m_matT.coeff(i+1,i+1)); 269*bf2c3715SXin Li RealScalar sd = numext::norm1(m_matT.coeff(i+1,i)); 270*bf2c3715SXin Li if (internal::isMuchSmallerThan(sd, d, NumTraits<RealScalar>::epsilon())) 271*bf2c3715SXin Li { 272*bf2c3715SXin Li m_matT.coeffRef(i+1,i) = ComplexScalar(0); 273*bf2c3715SXin Li return true; 274*bf2c3715SXin Li } 275*bf2c3715SXin Li return false; 276*bf2c3715SXin Li } 277*bf2c3715SXin Li 278*bf2c3715SXin Li 279*bf2c3715SXin Li /** Compute the shift in the current QR iteration. */ 280*bf2c3715SXin Li template<typename MatrixType> 281*bf2c3715SXin Li typename ComplexSchur<MatrixType>::ComplexScalar ComplexSchur<MatrixType>::computeShift(Index iu, Index iter) 282*bf2c3715SXin Li { 283*bf2c3715SXin Li using std::abs; 284*bf2c3715SXin Li if (iter == 10 || iter == 20) 285*bf2c3715SXin Li { 286*bf2c3715SXin Li // exceptional shift, taken from http://www.netlib.org/eispack/comqr.f 287*bf2c3715SXin Li return abs(numext::real(m_matT.coeff(iu,iu-1))) + abs(numext::real(m_matT.coeff(iu-1,iu-2))); 288*bf2c3715SXin Li } 289*bf2c3715SXin Li 290*bf2c3715SXin Li // compute the shift as one of the eigenvalues of t, the 2x2 291*bf2c3715SXin Li // diagonal block on the bottom of the active submatrix 292*bf2c3715SXin Li Matrix<ComplexScalar,2,2> t = m_matT.template block<2,2>(iu-1,iu-1); 293*bf2c3715SXin Li RealScalar normt = t.cwiseAbs().sum(); 294*bf2c3715SXin Li t /= normt; // the normalization by sf is to avoid under/overflow 295*bf2c3715SXin Li 296*bf2c3715SXin Li ComplexScalar b = t.coeff(0,1) * t.coeff(1,0); 297*bf2c3715SXin Li ComplexScalar c = t.coeff(0,0) - t.coeff(1,1); 298*bf2c3715SXin Li ComplexScalar disc = sqrt(c*c + RealScalar(4)*b); 299*bf2c3715SXin Li ComplexScalar det = t.coeff(0,0) * t.coeff(1,1) - b; 300*bf2c3715SXin Li ComplexScalar trace = t.coeff(0,0) + t.coeff(1,1); 301*bf2c3715SXin Li ComplexScalar eival1 = (trace + disc) / RealScalar(2); 302*bf2c3715SXin Li ComplexScalar eival2 = (trace - disc) / RealScalar(2); 303*bf2c3715SXin Li RealScalar eival1_norm = numext::norm1(eival1); 304*bf2c3715SXin Li RealScalar eival2_norm = numext::norm1(eival2); 305*bf2c3715SXin Li // A division by zero can only occur if eival1==eival2==0. 306*bf2c3715SXin Li // In this case, det==0, and all we have to do is checking that eival2_norm!=0 307*bf2c3715SXin Li if(eival1_norm > eival2_norm) 308*bf2c3715SXin Li eival2 = det / eival1; 309*bf2c3715SXin Li else if(eival2_norm!=RealScalar(0)) 310*bf2c3715SXin Li eival1 = det / eival2; 311*bf2c3715SXin Li 312*bf2c3715SXin Li // choose the eigenvalue closest to the bottom entry of the diagonal 313*bf2c3715SXin Li if(numext::norm1(eival1-t.coeff(1,1)) < numext::norm1(eival2-t.coeff(1,1))) 314*bf2c3715SXin Li return normt * eival1; 315*bf2c3715SXin Li else 316*bf2c3715SXin Li return normt * eival2; 317*bf2c3715SXin Li } 318*bf2c3715SXin Li 319*bf2c3715SXin Li 320*bf2c3715SXin Li template<typename MatrixType> 321*bf2c3715SXin Li template<typename InputType> 322*bf2c3715SXin Li ComplexSchur<MatrixType>& ComplexSchur<MatrixType>::compute(const EigenBase<InputType>& matrix, bool computeU) 323*bf2c3715SXin Li { 324*bf2c3715SXin Li m_matUisUptodate = false; 325*bf2c3715SXin Li eigen_assert(matrix.cols() == matrix.rows()); 326*bf2c3715SXin Li 327*bf2c3715SXin Li if(matrix.cols() == 1) 328*bf2c3715SXin Li { 329*bf2c3715SXin Li m_matT = matrix.derived().template cast<ComplexScalar>(); 330*bf2c3715SXin Li if(computeU) m_matU = ComplexMatrixType::Identity(1,1); 331*bf2c3715SXin Li m_info = Success; 332*bf2c3715SXin Li m_isInitialized = true; 333*bf2c3715SXin Li m_matUisUptodate = computeU; 334*bf2c3715SXin Li return *this; 335*bf2c3715SXin Li } 336*bf2c3715SXin Li 337*bf2c3715SXin Li internal::complex_schur_reduce_to_hessenberg<MatrixType, NumTraits<Scalar>::IsComplex>::run(*this, matrix.derived(), computeU); 338*bf2c3715SXin Li computeFromHessenberg(m_matT, m_matU, computeU); 339*bf2c3715SXin Li return *this; 340*bf2c3715SXin Li } 341*bf2c3715SXin Li 342*bf2c3715SXin Li template<typename MatrixType> 343*bf2c3715SXin Li template<typename HessMatrixType, typename OrthMatrixType> 344*bf2c3715SXin Li ComplexSchur<MatrixType>& ComplexSchur<MatrixType>::computeFromHessenberg(const HessMatrixType& matrixH, const OrthMatrixType& matrixQ, bool computeU) 345*bf2c3715SXin Li { 346*bf2c3715SXin Li m_matT = matrixH; 347*bf2c3715SXin Li if(computeU) 348*bf2c3715SXin Li m_matU = matrixQ; 349*bf2c3715SXin Li reduceToTriangularForm(computeU); 350*bf2c3715SXin Li return *this; 351*bf2c3715SXin Li } 352*bf2c3715SXin Li namespace internal { 353*bf2c3715SXin Li 354*bf2c3715SXin Li /* Reduce given matrix to Hessenberg form */ 355*bf2c3715SXin Li template<typename MatrixType, bool IsComplex> 356*bf2c3715SXin Li struct complex_schur_reduce_to_hessenberg 357*bf2c3715SXin Li { 358*bf2c3715SXin Li // this is the implementation for the case IsComplex = true 359*bf2c3715SXin Li static void run(ComplexSchur<MatrixType>& _this, const MatrixType& matrix, bool computeU) 360*bf2c3715SXin Li { 361*bf2c3715SXin Li _this.m_hess.compute(matrix); 362*bf2c3715SXin Li _this.m_matT = _this.m_hess.matrixH(); 363*bf2c3715SXin Li if(computeU) _this.m_matU = _this.m_hess.matrixQ(); 364*bf2c3715SXin Li } 365*bf2c3715SXin Li }; 366*bf2c3715SXin Li 367*bf2c3715SXin Li template<typename MatrixType> 368*bf2c3715SXin Li struct complex_schur_reduce_to_hessenberg<MatrixType, false> 369*bf2c3715SXin Li { 370*bf2c3715SXin Li static void run(ComplexSchur<MatrixType>& _this, const MatrixType& matrix, bool computeU) 371*bf2c3715SXin Li { 372*bf2c3715SXin Li typedef typename ComplexSchur<MatrixType>::ComplexScalar ComplexScalar; 373*bf2c3715SXin Li 374*bf2c3715SXin Li // Note: m_hess is over RealScalar; m_matT and m_matU is over ComplexScalar 375*bf2c3715SXin Li _this.m_hess.compute(matrix); 376*bf2c3715SXin Li _this.m_matT = _this.m_hess.matrixH().template cast<ComplexScalar>(); 377*bf2c3715SXin Li if(computeU) 378*bf2c3715SXin Li { 379*bf2c3715SXin Li // This may cause an allocation which seems to be avoidable 380*bf2c3715SXin Li MatrixType Q = _this.m_hess.matrixQ(); 381*bf2c3715SXin Li _this.m_matU = Q.template cast<ComplexScalar>(); 382*bf2c3715SXin Li } 383*bf2c3715SXin Li } 384*bf2c3715SXin Li }; 385*bf2c3715SXin Li 386*bf2c3715SXin Li } // end namespace internal 387*bf2c3715SXin Li 388*bf2c3715SXin Li // Reduce the Hessenberg matrix m_matT to triangular form by QR iteration. 389*bf2c3715SXin Li template<typename MatrixType> 390*bf2c3715SXin Li void ComplexSchur<MatrixType>::reduceToTriangularForm(bool computeU) 391*bf2c3715SXin Li { 392*bf2c3715SXin Li Index maxIters = m_maxIters; 393*bf2c3715SXin Li if (maxIters == -1) 394*bf2c3715SXin Li maxIters = m_maxIterationsPerRow * m_matT.rows(); 395*bf2c3715SXin Li 396*bf2c3715SXin Li // The matrix m_matT is divided in three parts. 397*bf2c3715SXin Li // Rows 0,...,il-1 are decoupled from the rest because m_matT(il,il-1) is zero. 398*bf2c3715SXin Li // Rows il,...,iu is the part we are working on (the active submatrix). 399*bf2c3715SXin Li // Rows iu+1,...,end are already brought in triangular form. 400*bf2c3715SXin Li Index iu = m_matT.cols() - 1; 401*bf2c3715SXin Li Index il; 402*bf2c3715SXin Li Index iter = 0; // number of iterations we are working on the (iu,iu) element 403*bf2c3715SXin Li Index totalIter = 0; // number of iterations for whole matrix 404*bf2c3715SXin Li 405*bf2c3715SXin Li while(true) 406*bf2c3715SXin Li { 407*bf2c3715SXin Li // find iu, the bottom row of the active submatrix 408*bf2c3715SXin Li while(iu > 0) 409*bf2c3715SXin Li { 410*bf2c3715SXin Li if(!subdiagonalEntryIsNeglegible(iu-1)) break; 411*bf2c3715SXin Li iter = 0; 412*bf2c3715SXin Li --iu; 413*bf2c3715SXin Li } 414*bf2c3715SXin Li 415*bf2c3715SXin Li // if iu is zero then we are done; the whole matrix is triangularized 416*bf2c3715SXin Li if(iu==0) break; 417*bf2c3715SXin Li 418*bf2c3715SXin Li // if we spent too many iterations, we give up 419*bf2c3715SXin Li iter++; 420*bf2c3715SXin Li totalIter++; 421*bf2c3715SXin Li if(totalIter > maxIters) break; 422*bf2c3715SXin Li 423*bf2c3715SXin Li // find il, the top row of the active submatrix 424*bf2c3715SXin Li il = iu-1; 425*bf2c3715SXin Li while(il > 0 && !subdiagonalEntryIsNeglegible(il-1)) 426*bf2c3715SXin Li { 427*bf2c3715SXin Li --il; 428*bf2c3715SXin Li } 429*bf2c3715SXin Li 430*bf2c3715SXin Li /* perform the QR step using Givens rotations. The first rotation 431*bf2c3715SXin Li creates a bulge; the (il+2,il) element becomes nonzero. This 432*bf2c3715SXin Li bulge is chased down to the bottom of the active submatrix. */ 433*bf2c3715SXin Li 434*bf2c3715SXin Li ComplexScalar shift = computeShift(iu, iter); 435*bf2c3715SXin Li JacobiRotation<ComplexScalar> rot; 436*bf2c3715SXin Li rot.makeGivens(m_matT.coeff(il,il) - shift, m_matT.coeff(il+1,il)); 437*bf2c3715SXin Li m_matT.rightCols(m_matT.cols()-il).applyOnTheLeft(il, il+1, rot.adjoint()); 438*bf2c3715SXin Li m_matT.topRows((std::min)(il+2,iu)+1).applyOnTheRight(il, il+1, rot); 439*bf2c3715SXin Li if(computeU) m_matU.applyOnTheRight(il, il+1, rot); 440*bf2c3715SXin Li 441*bf2c3715SXin Li for(Index i=il+1 ; i<iu ; i++) 442*bf2c3715SXin Li { 443*bf2c3715SXin Li rot.makeGivens(m_matT.coeffRef(i,i-1), m_matT.coeffRef(i+1,i-1), &m_matT.coeffRef(i,i-1)); 444*bf2c3715SXin Li m_matT.coeffRef(i+1,i-1) = ComplexScalar(0); 445*bf2c3715SXin Li m_matT.rightCols(m_matT.cols()-i).applyOnTheLeft(i, i+1, rot.adjoint()); 446*bf2c3715SXin Li m_matT.topRows((std::min)(i+2,iu)+1).applyOnTheRight(i, i+1, rot); 447*bf2c3715SXin Li if(computeU) m_matU.applyOnTheRight(i, i+1, rot); 448*bf2c3715SXin Li } 449*bf2c3715SXin Li } 450*bf2c3715SXin Li 451*bf2c3715SXin Li if(totalIter <= maxIters) 452*bf2c3715SXin Li m_info = Success; 453*bf2c3715SXin Li else 454*bf2c3715SXin Li m_info = NoConvergence; 455*bf2c3715SXin Li 456*bf2c3715SXin Li m_isInitialized = true; 457*bf2c3715SXin Li m_matUisUptodate = computeU; 458*bf2c3715SXin Li } 459*bf2c3715SXin Li 460*bf2c3715SXin Li } // end namespace Eigen 461*bf2c3715SXin Li 462*bf2c3715SXin Li #endif // EIGEN_COMPLEX_SCHUR_H 463