Important changes to repositories hosted on mbed.com
Mbed hosted mercurial repositories are deprecated and are due to be permanently deleted in July 2026.
To keep a copy of this software download the repository Zip archive or clone locally using Mercurial.
It is also possible to export all your personal repositories from the account settings page.
Tridiagonalization.h
00001 // This file is part of Eigen, a lightweight C++ template library 00002 // for linear algebra. 00003 // 00004 // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr> 00005 // Copyright (C) 2010 Jitse Niesen <jitse@maths.leeds.ac.uk> 00006 // 00007 // This Source Code Form is subject to the terms of the Mozilla 00008 // Public License v. 2.0. If a copy of the MPL was not distributed 00009 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. 00010 00011 #ifndef EIGEN_TRIDIAGONALIZATION_H 00012 #define EIGEN_TRIDIAGONALIZATION_H 00013 00014 namespace Eigen { 00015 00016 namespace internal { 00017 00018 template<typename MatrixType> struct TridiagonalizationMatrixTReturnType; 00019 template<typename MatrixType> 00020 struct traits<TridiagonalizationMatrixTReturnType<MatrixType> > 00021 { 00022 typedef typename MatrixType::PlainObject ReturnType; 00023 }; 00024 00025 template<typename MatrixType, typename CoeffVectorType> 00026 void tridiagonalization_inplace(MatrixType& matA, CoeffVectorType& hCoeffs); 00027 } 00028 00029 /** \eigenvalues_module \ingroup Eigenvalues_Module 00030 * 00031 * 00032 * \class Tridiagonalization 00033 * 00034 * \brief Tridiagonal decomposition of a selfadjoint matrix 00035 * 00036 * \tparam _MatrixType the type of the matrix of which we are computing the 00037 * tridiagonal decomposition; this is expected to be an instantiation of the 00038 * Matrix class template. 00039 * 00040 * This class performs a tridiagonal decomposition of a selfadjoint matrix \f$ A \f$ such that: 00041 * \f$ A = Q T Q^* \f$ where \f$ Q \f$ is unitary and \f$ T \f$ a real symmetric tridiagonal matrix. 00042 * 00043 * A tridiagonal matrix is a matrix which has nonzero elements only on the 00044 * main diagonal and the first diagonal below and above it. The Hessenberg 00045 * decomposition of a selfadjoint matrix is in fact a tridiagonal 00046 * decomposition. This class is used in SelfAdjointEigenSolver to compute the 00047 * eigenvalues and eigenvectors of a selfadjoint matrix. 00048 * 00049 * Call the function compute() to compute the tridiagonal decomposition of a 00050 * given matrix. Alternatively, you can use the Tridiagonalization(const MatrixType&) 00051 * constructor which computes the tridiagonal Schur decomposition at 00052 * construction time. Once the decomposition is computed, you can use the 00053 * matrixQ() and matrixT() functions to retrieve the matrices Q and T in the 00054 * decomposition. 00055 * 00056 * The documentation of Tridiagonalization(const MatrixType&) contains an 00057 * example of the typical use of this class. 00058 * 00059 * \sa class HessenbergDecomposition, class SelfAdjointEigenSolver 00060 */ 00061 template<typename _MatrixType> class Tridiagonalization 00062 { 00063 public: 00064 00065 /** \brief Synonym for the template parameter \p _MatrixType. */ 00066 typedef _MatrixType MatrixType; 00067 00068 typedef typename MatrixType::Scalar Scalar; 00069 typedef typename NumTraits<Scalar>::Real RealScalar; 00070 typedef typename MatrixType::Index Index; 00071 00072 enum { 00073 Size = MatrixType::RowsAtCompileTime, 00074 SizeMinusOne = Size == Dynamic ? Dynamic : (Size > 1 ? Size - 1 : 1), 00075 Options = MatrixType::Options, 00076 MaxSize = MatrixType::MaxRowsAtCompileTime, 00077 MaxSizeMinusOne = MaxSize == Dynamic ? Dynamic : (MaxSize > 1 ? MaxSize - 1 : 1) 00078 }; 00079 00080 typedef Matrix<Scalar, SizeMinusOne, 1, Options & ~RowMajor, MaxSizeMinusOne, 1> CoeffVectorType; 00081 typedef typename internal::plain_col_type<MatrixType, RealScalar>::type DiagonalType; 00082 typedef Matrix<RealScalar, SizeMinusOne, 1, Options & ~RowMajor, MaxSizeMinusOne, 1> SubDiagonalType; 00083 typedef typename internal::remove_all<typename MatrixType::RealReturnType>::type MatrixTypeRealView; 00084 typedef internal::TridiagonalizationMatrixTReturnType<MatrixTypeRealView> MatrixTReturnType; 00085 00086 typedef typename internal::conditional<NumTraits<Scalar>::IsComplex, 00087 typename internal::add_const_on_value_type<typename Diagonal<const MatrixType>::RealReturnType>::type, 00088 const Diagonal<const MatrixType> 00089 >::type DiagonalReturnType; 00090 00091 typedef typename internal::conditional<NumTraits<Scalar>::IsComplex, 00092 typename internal::add_const_on_value_type<typename Diagonal< 00093 Block<const MatrixType,SizeMinusOne,SizeMinusOne> >::RealReturnType>::type, 00094 const Diagonal< 00095 Block<const MatrixType,SizeMinusOne,SizeMinusOne> > 00096 >::type SubDiagonalReturnType; 00097 00098 /** \brief Return type of matrixQ() */ 00099 typedef HouseholderSequence<MatrixType,typename internal::remove_all<typename CoeffVectorType::ConjugateReturnType>::type > HouseholderSequenceType; 00100 00101 /** \brief Default constructor. 00102 * 00103 * \param [in] size Positive integer, size of the matrix whose tridiagonal 00104 * decomposition will be computed. 00105 * 00106 * The default constructor is useful in cases in which the user intends to 00107 * perform decompositions via compute(). The \p size parameter is only 00108 * used as a hint. It is not an error to give a wrong \p size, but it may 00109 * impair performance. 00110 * 00111 * \sa compute() for an example. 00112 */ 00113 Tridiagonalization(Index size = Size==Dynamic ? 2 : Size) 00114 : m_matrix(size,size), 00115 m_hCoeffs(size > 1 ? size-1 : 1), 00116 m_isInitialized(false) 00117 {} 00118 00119 /** \brief Constructor; computes tridiagonal decomposition of given matrix. 00120 * 00121 * \param[in] matrix Selfadjoint matrix whose tridiagonal decomposition 00122 * is to be computed. 00123 * 00124 * This constructor calls compute() to compute the tridiagonal decomposition. 00125 * 00126 * Example: \include Tridiagonalization_Tridiagonalization_MatrixType.cpp 00127 * Output: \verbinclude Tridiagonalization_Tridiagonalization_MatrixType.out 00128 */ 00129 Tridiagonalization(const MatrixType& matrix) 00130 : m_matrix(matrix), 00131 m_hCoeffs(matrix.cols() > 1 ? matrix.cols()-1 : 1), 00132 m_isInitialized(false) 00133 { 00134 internal::tridiagonalization_inplace(m_matrix, m_hCoeffs); 00135 m_isInitialized = true; 00136 } 00137 00138 /** \brief Computes tridiagonal decomposition of given matrix. 00139 * 00140 * \param[in] matrix Selfadjoint matrix whose tridiagonal decomposition 00141 * is to be computed. 00142 * \returns Reference to \c *this 00143 * 00144 * The tridiagonal decomposition is computed by bringing the columns of 00145 * the matrix successively in the required form using Householder 00146 * reflections. The cost is \f$ 4n^3/3 \f$ flops, where \f$ n \f$ denotes 00147 * the size of the given matrix. 00148 * 00149 * This method reuses of the allocated data in the Tridiagonalization 00150 * object, if the size of the matrix does not change. 00151 * 00152 * Example: \include Tridiagonalization_compute.cpp 00153 * Output: \verbinclude Tridiagonalization_compute.out 00154 */ 00155 Tridiagonalization & compute(const MatrixType& matrix) 00156 { 00157 m_matrix = matrix; 00158 m_hCoeffs.resize(matrix.rows()-1, 1); 00159 internal::tridiagonalization_inplace(m_matrix, m_hCoeffs); 00160 m_isInitialized = true; 00161 return *this; 00162 } 00163 00164 /** \brief Returns the Householder coefficients. 00165 * 00166 * \returns a const reference to the vector of Householder coefficients 00167 * 00168 * \pre Either the constructor Tridiagonalization(const MatrixType&) or 00169 * the member function compute(const MatrixType&) has been called before 00170 * to compute the tridiagonal decomposition of a matrix. 00171 * 00172 * The Householder coefficients allow the reconstruction of the matrix 00173 * \f$ Q \f$ in the tridiagonal decomposition from the packed data. 00174 * 00175 * Example: \include Tridiagonalization_householderCoefficients.cpp 00176 * Output: \verbinclude Tridiagonalization_householderCoefficients.out 00177 * 00178 * \sa packedMatrix(), \ref Householder_Module "Householder module" 00179 */ 00180 inline CoeffVectorType householderCoefficients() const 00181 { 00182 eigen_assert(m_isInitialized && "Tridiagonalization is not initialized."); 00183 return m_hCoeffs; 00184 } 00185 00186 /** \brief Returns the internal representation of the decomposition 00187 * 00188 * \returns a const reference to a matrix with the internal representation 00189 * of the decomposition. 00190 * 00191 * \pre Either the constructor Tridiagonalization(const MatrixType&) or 00192 * the member function compute(const MatrixType&) has been called before 00193 * to compute the tridiagonal decomposition of a matrix. 00194 * 00195 * The returned matrix contains the following information: 00196 * - the strict upper triangular part is equal to the input matrix A. 00197 * - the diagonal and lower sub-diagonal represent the real tridiagonal 00198 * symmetric matrix T. 00199 * - the rest of the lower part contains the Householder vectors that, 00200 * combined with Householder coefficients returned by 00201 * householderCoefficients(), allows to reconstruct the matrix Q as 00202 * \f$ Q = H_{N-1} \ldots H_1 H_0 \f$. 00203 * Here, the matrices \f$ H_i \f$ are the Householder transformations 00204 * \f$ H_i = (I - h_i v_i v_i^T) \f$ 00205 * where \f$ h_i \f$ is the \f$ i \f$th Householder coefficient and 00206 * \f$ v_i \f$ is the Householder vector defined by 00207 * \f$ v_i = [ 0, \ldots, 0, 1, M(i+2,i), \ldots, M(N-1,i) ]^T \f$ 00208 * with M the matrix returned by this function. 00209 * 00210 * See LAPACK for further details on this packed storage. 00211 * 00212 * Example: \include Tridiagonalization_packedMatrix.cpp 00213 * Output: \verbinclude Tridiagonalization_packedMatrix.out 00214 * 00215 * \sa householderCoefficients() 00216 */ 00217 inline const MatrixType& packedMatrix() const 00218 { 00219 eigen_assert(m_isInitialized && "Tridiagonalization is not initialized."); 00220 return m_matrix; 00221 } 00222 00223 /** \brief Returns the unitary matrix Q in the decomposition 00224 * 00225 * \returns object representing the matrix Q 00226 * 00227 * \pre Either the constructor Tridiagonalization(const MatrixType&) or 00228 * the member function compute(const MatrixType&) has been called before 00229 * to compute the tridiagonal decomposition of a matrix. 00230 * 00231 * This function returns a light-weight object of template class 00232 * HouseholderSequence. You can either apply it directly to a matrix or 00233 * you can convert it to a matrix of type #MatrixType. 00234 * 00235 * \sa Tridiagonalization(const MatrixType&) for an example, 00236 * matrixT(), class HouseholderSequence 00237 */ 00238 HouseholderSequenceType matrixQ() const 00239 { 00240 eigen_assert(m_isInitialized && "Tridiagonalization is not initialized."); 00241 return HouseholderSequenceType(m_matrix, m_hCoeffs.conjugate()) 00242 .setLength(m_matrix.rows() - 1) 00243 .setShift(1); 00244 } 00245 00246 /** \brief Returns an expression of the tridiagonal matrix T in the decomposition 00247 * 00248 * \returns expression object representing the matrix T 00249 * 00250 * \pre Either the constructor Tridiagonalization(const MatrixType&) or 00251 * the member function compute(const MatrixType&) has been called before 00252 * to compute the tridiagonal decomposition of a matrix. 00253 * 00254 * Currently, this function can be used to extract the matrix T from internal 00255 * data and copy it to a dense matrix object. In most cases, it may be 00256 * sufficient to directly use the packed matrix or the vector expressions 00257 * returned by diagonal() and subDiagonal() instead of creating a new 00258 * dense copy matrix with this function. 00259 * 00260 * \sa Tridiagonalization(const MatrixType&) for an example, 00261 * matrixQ(), packedMatrix(), diagonal(), subDiagonal() 00262 */ 00263 MatrixTReturnType matrixT() const 00264 { 00265 eigen_assert(m_isInitialized && "Tridiagonalization is not initialized."); 00266 return MatrixTReturnType(m_matrix.real()); 00267 } 00268 00269 /** \brief Returns the diagonal of the tridiagonal matrix T in the decomposition. 00270 * 00271 * \returns expression representing the diagonal of T 00272 * 00273 * \pre Either the constructor Tridiagonalization(const MatrixType&) or 00274 * the member function compute(const MatrixType&) has been called before 00275 * to compute the tridiagonal decomposition of a matrix. 00276 * 00277 * Example: \include Tridiagonalization_diagonal.cpp 00278 * Output: \verbinclude Tridiagonalization_diagonal.out 00279 * 00280 * \sa matrixT(), subDiagonal() 00281 */ 00282 DiagonalReturnType diagonal() const; 00283 00284 /** \brief Returns the subdiagonal of the tridiagonal matrix T in the decomposition. 00285 * 00286 * \returns expression representing the subdiagonal of T 00287 * 00288 * \pre Either the constructor Tridiagonalization(const MatrixType&) or 00289 * the member function compute(const MatrixType&) has been called before 00290 * to compute the tridiagonal decomposition of a matrix. 00291 * 00292 * \sa diagonal() for an example, matrixT() 00293 */ 00294 SubDiagonalReturnType subDiagonal() const; 00295 00296 protected: 00297 00298 MatrixType m_matrix; 00299 CoeffVectorType m_hCoeffs; 00300 bool m_isInitialized; 00301 }; 00302 00303 template<typename MatrixType> 00304 typename Tridiagonalization<MatrixType>::DiagonalReturnType 00305 Tridiagonalization<MatrixType>::diagonal() const 00306 { 00307 eigen_assert(m_isInitialized && "Tridiagonalization is not initialized."); 00308 return m_matrix.diagonal(); 00309 } 00310 00311 template<typename MatrixType> 00312 typename Tridiagonalization<MatrixType>::SubDiagonalReturnType 00313 Tridiagonalization<MatrixType>::subDiagonal() const 00314 { 00315 eigen_assert(m_isInitialized && "Tridiagonalization is not initialized."); 00316 Index n = m_matrix.rows(); 00317 return Block<const MatrixType,SizeMinusOne,SizeMinusOne>(m_matrix, 1, 0, n-1,n-1).diagonal(); 00318 } 00319 00320 namespace internal { 00321 00322 /** \internal 00323 * Performs a tridiagonal decomposition of the selfadjoint matrix \a matA in-place. 00324 * 00325 * \param[in,out] matA On input the selfadjoint matrix. Only the \b lower triangular part is referenced. 00326 * On output, the strict upper part is left unchanged, and the lower triangular part 00327 * represents the T and Q matrices in packed format has detailed below. 00328 * \param[out] hCoeffs returned Householder coefficients (see below) 00329 * 00330 * On output, the tridiagonal selfadjoint matrix T is stored in the diagonal 00331 * and lower sub-diagonal of the matrix \a matA. 00332 * The unitary matrix Q is represented in a compact way as a product of 00333 * Householder reflectors \f$ H_i \f$ such that: 00334 * \f$ Q = H_{N-1} \ldots H_1 H_0 \f$. 00335 * The Householder reflectors are defined as 00336 * \f$ H_i = (I - h_i v_i v_i^T) \f$ 00337 * where \f$ h_i = hCoeffs[i]\f$ is the \f$ i \f$th Householder coefficient and 00338 * \f$ v_i \f$ is the Householder vector defined by 00339 * \f$ v_i = [ 0, \ldots, 0, 1, matA(i+2,i), \ldots, matA(N-1,i) ]^T \f$. 00340 * 00341 * Implemented from Golub's "Matrix Computations", algorithm 8.3.1. 00342 * 00343 * \sa Tridiagonalization::packedMatrix() 00344 */ 00345 template<typename MatrixType, typename CoeffVectorType> 00346 void tridiagonalization_inplace(MatrixType& matA, CoeffVectorType& hCoeffs) 00347 { 00348 using numext::conj; 00349 typedef typename MatrixType::Index Index; 00350 typedef typename MatrixType::Scalar Scalar; 00351 typedef typename MatrixType::RealScalar RealScalar; 00352 Index n = matA.rows(); 00353 eigen_assert(n==matA.cols()); 00354 eigen_assert(n==hCoeffs.size()+1 || n==1); 00355 00356 for (Index i = 0; i<n-1; ++i) 00357 { 00358 Index remainingSize = n-i-1; 00359 RealScalar beta; 00360 Scalar h; 00361 matA.col(i).tail(remainingSize).makeHouseholderInPlace(h, beta); 00362 00363 // Apply similarity transformation to remaining columns, 00364 // i.e., A = H A H' where H = I - h v v' and v = matA.col(i).tail(n-i-1) 00365 matA.col(i).coeffRef(i+1) = 1; 00366 00367 hCoeffs.tail(n-i-1).noalias() = (matA.bottomRightCorner(remainingSize,remainingSize).template selfadjointView<Lower>() 00368 * (conj(h) * matA.col(i).tail(remainingSize))); 00369 00370 hCoeffs.tail(n-i-1) += (conj(h)*Scalar(-0.5)*(hCoeffs.tail(remainingSize).dot(matA.col(i).tail(remainingSize)))) * matA.col(i).tail(n-i-1); 00371 00372 matA.bottomRightCorner(remainingSize, remainingSize).template selfadjointView<Lower>() 00373 .rankUpdate(matA.col(i).tail(remainingSize), hCoeffs.tail(remainingSize), -1); 00374 00375 matA.col(i).coeffRef(i+1) = beta; 00376 hCoeffs.coeffRef(i) = h; 00377 } 00378 } 00379 00380 // forward declaration, implementation at the end of this file 00381 template<typename MatrixType, 00382 int Size=MatrixType::ColsAtCompileTime, 00383 bool IsComplex=NumTraits<typename MatrixType::Scalar>::IsComplex> 00384 struct tridiagonalization_inplace_selector; 00385 00386 /** \brief Performs a full tridiagonalization in place 00387 * 00388 * \param[in,out] mat On input, the selfadjoint matrix whose tridiagonal 00389 * decomposition is to be computed. Only the lower triangular part referenced. 00390 * The rest is left unchanged. On output, the orthogonal matrix Q 00391 * in the decomposition if \p extractQ is true. 00392 * \param[out] diag The diagonal of the tridiagonal matrix T in the 00393 * decomposition. 00394 * \param[out] subdiag The subdiagonal of the tridiagonal matrix T in 00395 * the decomposition. 00396 * \param[in] extractQ If true, the orthogonal matrix Q in the 00397 * decomposition is computed and stored in \p mat. 00398 * 00399 * Computes the tridiagonal decomposition of the selfadjoint matrix \p mat in place 00400 * such that \f$ mat = Q T Q^* \f$ where \f$ Q \f$ is unitary and \f$ T \f$ a real 00401 * symmetric tridiagonal matrix. 00402 * 00403 * The tridiagonal matrix T is passed to the output parameters \p diag and \p subdiag. If 00404 * \p extractQ is true, then the orthogonal matrix Q is passed to \p mat. Otherwise the lower 00405 * part of the matrix \p mat is destroyed. 00406 * 00407 * The vectors \p diag and \p subdiag are not resized. The function 00408 * assumes that they are already of the correct size. The length of the 00409 * vector \p diag should equal the number of rows in \p mat, and the 00410 * length of the vector \p subdiag should be one left. 00411 * 00412 * This implementation contains an optimized path for 3-by-3 matrices 00413 * which is especially useful for plane fitting. 00414 * 00415 * \note Currently, it requires two temporary vectors to hold the intermediate 00416 * Householder coefficients, and to reconstruct the matrix Q from the Householder 00417 * reflectors. 00418 * 00419 * Example (this uses the same matrix as the example in 00420 * Tridiagonalization::Tridiagonalization(const MatrixType&)): 00421 * \include Tridiagonalization_decomposeInPlace.cpp 00422 * Output: \verbinclude Tridiagonalization_decomposeInPlace.out 00423 * 00424 * \sa class Tridiagonalization 00425 */ 00426 template<typename MatrixType, typename DiagonalType, typename SubDiagonalType> 00427 void tridiagonalization_inplace(MatrixType& mat, DiagonalType& diag, SubDiagonalType& subdiag, bool extractQ) 00428 { 00429 eigen_assert(mat.cols()==mat.rows() && diag.size()==mat.rows() && subdiag.size()==mat.rows()-1); 00430 tridiagonalization_inplace_selector<MatrixType>::run(mat, diag, subdiag, extractQ); 00431 } 00432 00433 /** \internal 00434 * General full tridiagonalization 00435 */ 00436 template<typename MatrixType, int Size, bool IsComplex> 00437 struct tridiagonalization_inplace_selector 00438 { 00439 typedef typename Tridiagonalization<MatrixType>::CoeffVectorType CoeffVectorType; 00440 typedef typename Tridiagonalization<MatrixType>::HouseholderSequenceType HouseholderSequenceType; 00441 typedef typename MatrixType::Index Index; 00442 template<typename DiagonalType, typename SubDiagonalType> 00443 static void run(MatrixType& mat, DiagonalType& diag, SubDiagonalType& subdiag, bool extractQ) 00444 { 00445 CoeffVectorType hCoeffs(mat.cols()-1); 00446 tridiagonalization_inplace(mat,hCoeffs); 00447 diag = mat.diagonal().real(); 00448 subdiag = mat.template diagonal<-1>().real(); 00449 if(extractQ) 00450 mat = HouseholderSequenceType(mat, hCoeffs.conjugate()) 00451 .setLength(mat.rows() - 1) 00452 .setShift(1); 00453 } 00454 }; 00455 00456 /** \internal 00457 * Specialization for 3x3 real matrices. 00458 * Especially useful for plane fitting. 00459 */ 00460 template<typename MatrixType> 00461 struct tridiagonalization_inplace_selector<MatrixType,3,false> 00462 { 00463 typedef typename MatrixType::Scalar Scalar; 00464 typedef typename MatrixType::RealScalar RealScalar; 00465 00466 template<typename DiagonalType, typename SubDiagonalType> 00467 static void run(MatrixType& mat, DiagonalType& diag, SubDiagonalType& subdiag, bool extractQ) 00468 { 00469 using std::sqrt; 00470 diag[0] = mat(0,0); 00471 RealScalar v1norm2 = numext::abs2(mat(2,0)); 00472 if(v1norm2 == RealScalar(0)) 00473 { 00474 diag[1] = mat(1,1); 00475 diag[2] = mat(2,2); 00476 subdiag[0] = mat(1,0); 00477 subdiag[1] = mat(2,1); 00478 if (extractQ) 00479 mat.setIdentity(); 00480 } 00481 else 00482 { 00483 RealScalar beta = sqrt(numext::abs2(mat(1,0)) + v1norm2); 00484 RealScalar invBeta = RealScalar(1)/beta; 00485 Scalar m01 = mat(1,0) * invBeta; 00486 Scalar m02 = mat(2,0) * invBeta; 00487 Scalar q = RealScalar(2)*m01*mat(2,1) + m02*(mat(2,2) - mat(1,1)); 00488 diag[1] = mat(1,1) + m02*q; 00489 diag[2] = mat(2,2) - m02*q; 00490 subdiag[0] = beta; 00491 subdiag[1] = mat(2,1) - m01 * q; 00492 if (extractQ) 00493 { 00494 mat << 1, 0, 0, 00495 0, m01, m02, 00496 0, m02, -m01; 00497 } 00498 } 00499 } 00500 }; 00501 00502 /** \internal 00503 * Trivial specialization for 1x1 matrices 00504 */ 00505 template<typename MatrixType, bool IsComplex> 00506 struct tridiagonalization_inplace_selector<MatrixType,1,IsComplex> 00507 { 00508 typedef typename MatrixType::Scalar Scalar; 00509 00510 template<typename DiagonalType, typename SubDiagonalType> 00511 static void run(MatrixType& mat, DiagonalType& diag, SubDiagonalType&, bool extractQ) 00512 { 00513 diag(0,0) = numext::real(mat(0,0)); 00514 if(extractQ) 00515 mat(0,0) = Scalar(1); 00516 } 00517 }; 00518 00519 /** \internal 00520 * \eigenvalues_module \ingroup Eigenvalues_Module 00521 * 00522 * \brief Expression type for return value of Tridiagonalization::matrixT() 00523 * 00524 * \tparam MatrixType type of underlying dense matrix 00525 */ 00526 template<typename MatrixType> struct TridiagonalizationMatrixTReturnType 00527 : public ReturnByValue<TridiagonalizationMatrixTReturnType<MatrixType> > 00528 { 00529 typedef typename MatrixType::Index Index; 00530 public: 00531 /** \brief Constructor. 00532 * 00533 * \param[in] mat The underlying dense matrix 00534 */ 00535 TridiagonalizationMatrixTReturnType(const MatrixType& mat) : m_matrix(mat) { } 00536 00537 template <typename ResultType> 00538 inline void evalTo(ResultType& result) const 00539 { 00540 result.setZero(); 00541 result.template diagonal<1>() = m_matrix.template diagonal<-1>().conjugate(); 00542 result.diagonal() = m_matrix.diagonal(); 00543 result.template diagonal<-1>() = m_matrix.template diagonal<-1>(); 00544 } 00545 00546 Index rows() const { return m_matrix.rows(); } 00547 Index cols() const { return m_matrix.cols(); } 00548 00549 protected: 00550 typename MatrixType::Nested m_matrix; 00551 }; 00552 00553 } // end namespace internal 00554 00555 } // end namespace Eigen 00556 00557 #endif // EIGEN_TRIDIAGONALIZATION_H
Generated on Thu Nov 17 2022 22:01:30 by
1.7.2