• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
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 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 
11 #ifndef EIGEN_TRIDIAGONALIZATION_H
12 #define EIGEN_TRIDIAGONALIZATION_H
13 
14 namespace Eigen {
15 
16 namespace internal {
17 
18 template<typename MatrixType> struct TridiagonalizationMatrixTReturnType;
19 template<typename MatrixType>
20 struct traits<TridiagonalizationMatrixTReturnType<MatrixType> >
21   : public traits<typename MatrixType::PlainObject>
22 {
23   typedef typename MatrixType::PlainObject ReturnType; // FIXME shall it be a BandMatrix?
24   enum { Flags = 0 };
25 };
26 
27 template<typename MatrixType, typename CoeffVectorType>
28 void tridiagonalization_inplace(MatrixType& matA, CoeffVectorType& hCoeffs);
29 }
30 
31 /** \eigenvalues_module \ingroup Eigenvalues_Module
32   *
33   *
34   * \class Tridiagonalization
35   *
36   * \brief Tridiagonal decomposition of a selfadjoint matrix
37   *
38   * \tparam _MatrixType the type of the matrix of which we are computing the
39   * tridiagonal decomposition; this is expected to be an instantiation of the
40   * Matrix class template.
41   *
42   * This class performs a tridiagonal decomposition of a selfadjoint matrix \f$ A \f$ such that:
43   * \f$ A = Q T Q^* \f$ where \f$ Q \f$ is unitary and \f$ T \f$ a real symmetric tridiagonal matrix.
44   *
45   * A tridiagonal matrix is a matrix which has nonzero elements only on the
46   * main diagonal and the first diagonal below and above it. The Hessenberg
47   * decomposition of a selfadjoint matrix is in fact a tridiagonal
48   * decomposition. This class is used in SelfAdjointEigenSolver to compute the
49   * eigenvalues and eigenvectors of a selfadjoint matrix.
50   *
51   * Call the function compute() to compute the tridiagonal decomposition of a
52   * given matrix. Alternatively, you can use the Tridiagonalization(const MatrixType&)
53   * constructor which computes the tridiagonal Schur decomposition at
54   * construction time. Once the decomposition is computed, you can use the
55   * matrixQ() and matrixT() functions to retrieve the matrices Q and T in the
56   * decomposition.
57   *
58   * The documentation of Tridiagonalization(const MatrixType&) contains an
59   * example of the typical use of this class.
60   *
61   * \sa class HessenbergDecomposition, class SelfAdjointEigenSolver
62   */
63 template<typename _MatrixType> class Tridiagonalization
64 {
65   public:
66 
67     /** \brief Synonym for the template parameter \p _MatrixType. */
68     typedef _MatrixType MatrixType;
69 
70     typedef typename MatrixType::Scalar Scalar;
71     typedef typename NumTraits<Scalar>::Real RealScalar;
72     typedef Eigen::Index Index; ///< \deprecated since Eigen 3.3
73 
74     enum {
75       Size = MatrixType::RowsAtCompileTime,
76       SizeMinusOne = Size == Dynamic ? Dynamic : (Size > 1 ? Size - 1 : 1),
77       Options = MatrixType::Options,
78       MaxSize = MatrixType::MaxRowsAtCompileTime,
79       MaxSizeMinusOne = MaxSize == Dynamic ? Dynamic : (MaxSize > 1 ? MaxSize - 1 : 1)
80     };
81 
82     typedef Matrix<Scalar, SizeMinusOne, 1, Options & ~RowMajor, MaxSizeMinusOne, 1> CoeffVectorType;
83     typedef typename internal::plain_col_type<MatrixType, RealScalar>::type DiagonalType;
84     typedef Matrix<RealScalar, SizeMinusOne, 1, Options & ~RowMajor, MaxSizeMinusOne, 1> SubDiagonalType;
85     typedef typename internal::remove_all<typename MatrixType::RealReturnType>::type MatrixTypeRealView;
86     typedef internal::TridiagonalizationMatrixTReturnType<MatrixTypeRealView> MatrixTReturnType;
87 
88     typedef typename internal::conditional<NumTraits<Scalar>::IsComplex,
89               typename internal::add_const_on_value_type<typename Diagonal<const MatrixType>::RealReturnType>::type,
90               const Diagonal<const MatrixType>
91             >::type DiagonalReturnType;
92 
93     typedef typename internal::conditional<NumTraits<Scalar>::IsComplex,
94               typename internal::add_const_on_value_type<typename Diagonal<const MatrixType, -1>::RealReturnType>::type,
95               const Diagonal<const MatrixType, -1>
96             >::type SubDiagonalReturnType;
97 
98     /** \brief Return type of matrixQ() */
99     typedef HouseholderSequence<MatrixType,typename internal::remove_all<typename CoeffVectorType::ConjugateReturnType>::type> HouseholderSequenceType;
100 
101     /** \brief Default constructor.
102       *
103       * \param [in]  size  Positive integer, size of the matrix whose tridiagonal
104       * decomposition will be computed.
105       *
106       * The default constructor is useful in cases in which the user intends to
107       * perform decompositions via compute().  The \p size parameter is only
108       * used as a hint. It is not an error to give a wrong \p size, but it may
109       * impair performance.
110       *
111       * \sa compute() for an example.
112       */
113     explicit Tridiagonalization(Index size = Size==Dynamic ? 2 : Size)
114       : m_matrix(size,size),
115         m_hCoeffs(size > 1 ? size-1 : 1),
116         m_isInitialized(false)
117     {}
118 
119     /** \brief Constructor; computes tridiagonal decomposition of given matrix.
120       *
121       * \param[in]  matrix  Selfadjoint matrix whose tridiagonal decomposition
122       * is to be computed.
123       *
124       * This constructor calls compute() to compute the tridiagonal decomposition.
125       *
126       * Example: \include Tridiagonalization_Tridiagonalization_MatrixType.cpp
127       * Output: \verbinclude Tridiagonalization_Tridiagonalization_MatrixType.out
128       */
129     template<typename InputType>
130     explicit Tridiagonalization(const EigenBase<InputType>& matrix)
131       : m_matrix(matrix.derived()),
132         m_hCoeffs(matrix.cols() > 1 ? matrix.cols()-1 : 1),
133         m_isInitialized(false)
134     {
135       internal::tridiagonalization_inplace(m_matrix, m_hCoeffs);
136       m_isInitialized = true;
137     }
138 
139     /** \brief Computes tridiagonal decomposition of given matrix.
140       *
141       * \param[in]  matrix  Selfadjoint matrix whose tridiagonal decomposition
142       * is to be computed.
143       * \returns    Reference to \c *this
144       *
145       * The tridiagonal decomposition is computed by bringing the columns of
146       * the matrix successively in the required form using Householder
147       * reflections. The cost is \f$ 4n^3/3 \f$ flops, where \f$ n \f$ denotes
148       * the size of the given matrix.
149       *
150       * This method reuses of the allocated data in the Tridiagonalization
151       * object, if the size of the matrix does not change.
152       *
153       * Example: \include Tridiagonalization_compute.cpp
154       * Output: \verbinclude Tridiagonalization_compute.out
155       */
156     template<typename InputType>
157     Tridiagonalization& compute(const EigenBase<InputType>& matrix)
158     {
159       m_matrix = matrix.derived();
160       m_hCoeffs.resize(matrix.rows()-1, 1);
161       internal::tridiagonalization_inplace(m_matrix, m_hCoeffs);
162       m_isInitialized = true;
163       return *this;
164     }
165 
166     /** \brief Returns the Householder coefficients.
167       *
168       * \returns a const reference to the vector of Householder coefficients
169       *
170       * \pre Either the constructor Tridiagonalization(const MatrixType&) or
171       * the member function compute(const MatrixType&) has been called before
172       * to compute the tridiagonal decomposition of a matrix.
173       *
174       * The Householder coefficients allow the reconstruction of the matrix
175       * \f$ Q \f$ in the tridiagonal decomposition from the packed data.
176       *
177       * Example: \include Tridiagonalization_householderCoefficients.cpp
178       * Output: \verbinclude Tridiagonalization_householderCoefficients.out
179       *
180       * \sa packedMatrix(), \ref Householder_Module "Householder module"
181       */
182     inline CoeffVectorType householderCoefficients() const
183     {
184       eigen_assert(m_isInitialized && "Tridiagonalization is not initialized.");
185       return m_hCoeffs;
186     }
187 
188     /** \brief Returns the internal representation of the decomposition
189       *
190       *	\returns a const reference to a matrix with the internal representation
191       *	         of the decomposition.
192       *
193       * \pre Either the constructor Tridiagonalization(const MatrixType&) or
194       * the member function compute(const MatrixType&) has been called before
195       * to compute the tridiagonal decomposition of a matrix.
196       *
197       * The returned matrix contains the following information:
198       *  - the strict upper triangular part is equal to the input matrix A.
199       *  - the diagonal and lower sub-diagonal represent the real tridiagonal
200       *    symmetric matrix T.
201       *  - the rest of the lower part contains the Householder vectors that,
202       *    combined with Householder coefficients returned by
203       *    householderCoefficients(), allows to reconstruct the matrix Q as
204       *       \f$ Q = H_{N-1} \ldots H_1 H_0 \f$.
205       *    Here, the matrices \f$ H_i \f$ are the Householder transformations
206       *       \f$ H_i = (I - h_i v_i v_i^T) \f$
207       *    where \f$ h_i \f$ is the \f$ i \f$th Householder coefficient and
208       *    \f$ v_i \f$ is the Householder vector defined by
209       *       \f$ v_i = [ 0, \ldots, 0, 1, M(i+2,i), \ldots, M(N-1,i) ]^T \f$
210       *    with M the matrix returned by this function.
211       *
212       * See LAPACK for further details on this packed storage.
213       *
214       * Example: \include Tridiagonalization_packedMatrix.cpp
215       * Output: \verbinclude Tridiagonalization_packedMatrix.out
216       *
217       * \sa householderCoefficients()
218       */
219     inline const MatrixType& packedMatrix() const
220     {
221       eigen_assert(m_isInitialized && "Tridiagonalization is not initialized.");
222       return m_matrix;
223     }
224 
225     /** \brief Returns the unitary matrix Q in the decomposition
226       *
227       * \returns object representing the matrix Q
228       *
229       * \pre Either the constructor Tridiagonalization(const MatrixType&) or
230       * the member function compute(const MatrixType&) has been called before
231       * to compute the tridiagonal decomposition of a matrix.
232       *
233       * This function returns a light-weight object of template class
234       * HouseholderSequence. You can either apply it directly to a matrix or
235       * you can convert it to a matrix of type #MatrixType.
236       *
237       * \sa Tridiagonalization(const MatrixType&) for an example,
238       *     matrixT(), class HouseholderSequence
239       */
240     HouseholderSequenceType matrixQ() const
241     {
242       eigen_assert(m_isInitialized && "Tridiagonalization is not initialized.");
243       return HouseholderSequenceType(m_matrix, m_hCoeffs.conjugate())
244              .setLength(m_matrix.rows() - 1)
245              .setShift(1);
246     }
247 
248     /** \brief Returns an expression of the tridiagonal matrix T in the decomposition
249       *
250       * \returns expression object representing the matrix T
251       *
252       * \pre Either the constructor Tridiagonalization(const MatrixType&) or
253       * the member function compute(const MatrixType&) has been called before
254       * to compute the tridiagonal decomposition of a matrix.
255       *
256       * Currently, this function can be used to extract the matrix T from internal
257       * data and copy it to a dense matrix object. In most cases, it may be
258       * sufficient to directly use the packed matrix or the vector expressions
259       * returned by diagonal() and subDiagonal() instead of creating a new
260       * dense copy matrix with this function.
261       *
262       * \sa Tridiagonalization(const MatrixType&) for an example,
263       * matrixQ(), packedMatrix(), diagonal(), subDiagonal()
264       */
265     MatrixTReturnType matrixT() const
266     {
267       eigen_assert(m_isInitialized && "Tridiagonalization is not initialized.");
268       return MatrixTReturnType(m_matrix.real());
269     }
270 
271     /** \brief Returns the diagonal of the tridiagonal matrix T in the decomposition.
272       *
273       * \returns expression representing the diagonal of T
274       *
275       * \pre Either the constructor Tridiagonalization(const MatrixType&) or
276       * the member function compute(const MatrixType&) has been called before
277       * to compute the tridiagonal decomposition of a matrix.
278       *
279       * Example: \include Tridiagonalization_diagonal.cpp
280       * Output: \verbinclude Tridiagonalization_diagonal.out
281       *
282       * \sa matrixT(), subDiagonal()
283       */
284     DiagonalReturnType diagonal() const;
285 
286     /** \brief Returns the subdiagonal of the tridiagonal matrix T in the decomposition.
287       *
288       * \returns expression representing the subdiagonal of T
289       *
290       * \pre Either the constructor Tridiagonalization(const MatrixType&) or
291       * the member function compute(const MatrixType&) has been called before
292       * to compute the tridiagonal decomposition of a matrix.
293       *
294       * \sa diagonal() for an example, matrixT()
295       */
296     SubDiagonalReturnType subDiagonal() const;
297 
298   protected:
299 
300     MatrixType m_matrix;
301     CoeffVectorType m_hCoeffs;
302     bool m_isInitialized;
303 };
304 
305 template<typename MatrixType>
306 typename Tridiagonalization<MatrixType>::DiagonalReturnType
307 Tridiagonalization<MatrixType>::diagonal() const
308 {
309   eigen_assert(m_isInitialized && "Tridiagonalization is not initialized.");
310   return m_matrix.diagonal().real();
311 }
312 
313 template<typename MatrixType>
314 typename Tridiagonalization<MatrixType>::SubDiagonalReturnType
315 Tridiagonalization<MatrixType>::subDiagonal() const
316 {
317   eigen_assert(m_isInitialized && "Tridiagonalization is not initialized.");
318   return m_matrix.template diagonal<-1>().real();
319 }
320 
321 namespace internal {
322 
323 /** \internal
324   * Performs a tridiagonal decomposition of the selfadjoint matrix \a matA in-place.
325   *
326   * \param[in,out] matA On input the selfadjoint matrix. Only the \b lower triangular part is referenced.
327   *                     On output, the strict upper part is left unchanged, and the lower triangular part
328   *                     represents the T and Q matrices in packed format has detailed below.
329   * \param[out]    hCoeffs returned Householder coefficients (see below)
330   *
331   * On output, the tridiagonal selfadjoint matrix T is stored in the diagonal
332   * and lower sub-diagonal of the matrix \a matA.
333   * The unitary matrix Q is represented in a compact way as a product of
334   * Householder reflectors \f$ H_i \f$ such that:
335   *       \f$ Q = H_{N-1} \ldots H_1 H_0 \f$.
336   * The Householder reflectors are defined as
337   *       \f$ H_i = (I - h_i v_i v_i^T) \f$
338   * where \f$ h_i = hCoeffs[i]\f$ is the \f$ i \f$th Householder coefficient and
339   * \f$ v_i \f$ is the Householder vector defined by
340   *       \f$ v_i = [ 0, \ldots, 0, 1, matA(i+2,i), \ldots, matA(N-1,i) ]^T \f$.
341   *
342   * Implemented from Golub's "Matrix Computations", algorithm 8.3.1.
343   *
344   * \sa Tridiagonalization::packedMatrix()
345   */
346 template<typename MatrixType, typename CoeffVectorType>
347 void tridiagonalization_inplace(MatrixType& matA, CoeffVectorType& hCoeffs)
348 {
349   using numext::conj;
350   typedef typename MatrixType::Scalar Scalar;
351   typedef typename MatrixType::RealScalar RealScalar;
352   Index n = matA.rows();
353   eigen_assert(n==matA.cols());
354   eigen_assert(n==hCoeffs.size()+1 || n==1);
355 
356   for (Index i = 0; i<n-1; ++i)
357   {
358     Index remainingSize = n-i-1;
359     RealScalar beta;
360     Scalar h;
361     matA.col(i).tail(remainingSize).makeHouseholderInPlace(h, beta);
362 
363     // Apply similarity transformation to remaining columns,
364     // i.e., A = H A H' where H = I - h v v' and v = matA.col(i).tail(n-i-1)
365     matA.col(i).coeffRef(i+1) = 1;
366 
367     hCoeffs.tail(n-i-1).noalias() = (matA.bottomRightCorner(remainingSize,remainingSize).template selfadjointView<Lower>()
368                                   * (conj(h) * matA.col(i).tail(remainingSize)));
369 
370     hCoeffs.tail(n-i-1) += (conj(h)*RealScalar(-0.5)*(hCoeffs.tail(remainingSize).dot(matA.col(i).tail(remainingSize)))) * matA.col(i).tail(n-i-1);
371 
372     matA.bottomRightCorner(remainingSize, remainingSize).template selfadjointView<Lower>()
373       .rankUpdate(matA.col(i).tail(remainingSize), hCoeffs.tail(remainingSize), Scalar(-1));
374 
375     matA.col(i).coeffRef(i+1) = beta;
376     hCoeffs.coeffRef(i) = h;
377   }
378 }
379 
380 // forward declaration, implementation at the end of this file
381 template<typename MatrixType,
382          int Size=MatrixType::ColsAtCompileTime,
383          bool IsComplex=NumTraits<typename MatrixType::Scalar>::IsComplex>
384 struct tridiagonalization_inplace_selector;
385 
386 /** \brief Performs a full tridiagonalization in place
387   *
388   * \param[in,out]  mat  On input, the selfadjoint matrix whose tridiagonal
389   *    decomposition is to be computed. Only the lower triangular part referenced.
390   *    The rest is left unchanged. On output, the orthogonal matrix Q
391   *    in the decomposition if \p extractQ is true.
392   * \param[out]  diag  The diagonal of the tridiagonal matrix T in the
393   *    decomposition.
394   * \param[out]  subdiag  The subdiagonal of the tridiagonal matrix T in
395   *    the decomposition.
396   * \param[in]  extractQ  If true, the orthogonal matrix Q in the
397   *    decomposition is computed and stored in \p mat.
398   *
399   * Computes the tridiagonal decomposition of the selfadjoint matrix \p mat in place
400   * such that \f$ mat = Q T Q^* \f$ where \f$ Q \f$ is unitary and \f$ T \f$ a real
401   * symmetric tridiagonal matrix.
402   *
403   * The tridiagonal matrix T is passed to the output parameters \p diag and \p subdiag. If
404   * \p extractQ is true, then the orthogonal matrix Q is passed to \p mat. Otherwise the lower
405   * part of the matrix \p mat is destroyed.
406   *
407   * The vectors \p diag and \p subdiag are not resized. The function
408   * assumes that they are already of the correct size. The length of the
409   * vector \p diag should equal the number of rows in \p mat, and the
410   * length of the vector \p subdiag should be one left.
411   *
412   * This implementation contains an optimized path for 3-by-3 matrices
413   * which is especially useful for plane fitting.
414   *
415   * \note Currently, it requires two temporary vectors to hold the intermediate
416   * Householder coefficients, and to reconstruct the matrix Q from the Householder
417   * reflectors.
418   *
419   * Example (this uses the same matrix as the example in
420   *    Tridiagonalization::Tridiagonalization(const MatrixType&)):
421   *    \include Tridiagonalization_decomposeInPlace.cpp
422   * Output: \verbinclude Tridiagonalization_decomposeInPlace.out
423   *
424   * \sa class Tridiagonalization
425   */
426 template<typename MatrixType, typename DiagonalType, typename SubDiagonalType>
427 void tridiagonalization_inplace(MatrixType& mat, DiagonalType& diag, SubDiagonalType& subdiag, bool extractQ)
428 {
429   eigen_assert(mat.cols()==mat.rows() && diag.size()==mat.rows() && subdiag.size()==mat.rows()-1);
430   tridiagonalization_inplace_selector<MatrixType>::run(mat, diag, subdiag, extractQ);
431 }
432 
433 /** \internal
434   * General full tridiagonalization
435   */
436 template<typename MatrixType, int Size, bool IsComplex>
437 struct tridiagonalization_inplace_selector
438 {
439   typedef typename Tridiagonalization<MatrixType>::CoeffVectorType CoeffVectorType;
440   typedef typename Tridiagonalization<MatrixType>::HouseholderSequenceType HouseholderSequenceType;
441   template<typename DiagonalType, typename SubDiagonalType>
442   static void run(MatrixType& mat, DiagonalType& diag, SubDiagonalType& subdiag, bool extractQ)
443   {
444     CoeffVectorType hCoeffs(mat.cols()-1);
445     tridiagonalization_inplace(mat,hCoeffs);
446     diag = mat.diagonal().real();
447     subdiag = mat.template diagonal<-1>().real();
448     if(extractQ)
449       mat = HouseholderSequenceType(mat, hCoeffs.conjugate())
450             .setLength(mat.rows() - 1)
451             .setShift(1);
452   }
453 };
454 
455 /** \internal
456   * Specialization for 3x3 real matrices.
457   * Especially useful for plane fitting.
458   */
459 template<typename MatrixType>
460 struct tridiagonalization_inplace_selector<MatrixType,3,false>
461 {
462   typedef typename MatrixType::Scalar Scalar;
463   typedef typename MatrixType::RealScalar RealScalar;
464 
465   template<typename DiagonalType, typename SubDiagonalType>
466   static void run(MatrixType& mat, DiagonalType& diag, SubDiagonalType& subdiag, bool extractQ)
467   {
468     using std::sqrt;
469     const RealScalar tol = (std::numeric_limits<RealScalar>::min)();
470     diag[0] = mat(0,0);
471     RealScalar v1norm2 = numext::abs2(mat(2,0));
472     if(v1norm2 <= tol)
473     {
474       diag[1] = mat(1,1);
475       diag[2] = mat(2,2);
476       subdiag[0] = mat(1,0);
477       subdiag[1] = mat(2,1);
478       if (extractQ)
479         mat.setIdentity();
480     }
481     else
482     {
483       RealScalar beta = sqrt(numext::abs2(mat(1,0)) + v1norm2);
484       RealScalar invBeta = RealScalar(1)/beta;
485       Scalar m01 = mat(1,0) * invBeta;
486       Scalar m02 = mat(2,0) * invBeta;
487       Scalar q = RealScalar(2)*m01*mat(2,1) + m02*(mat(2,2) - mat(1,1));
488       diag[1] = mat(1,1) + m02*q;
489       diag[2] = mat(2,2) - m02*q;
490       subdiag[0] = beta;
491       subdiag[1] = mat(2,1) - m01 * q;
492       if (extractQ)
493       {
494         mat << 1,   0,    0,
495                0, m01,  m02,
496                0, m02, -m01;
497       }
498     }
499   }
500 };
501 
502 /** \internal
503   * Trivial specialization for 1x1 matrices
504   */
505 template<typename MatrixType, bool IsComplex>
506 struct tridiagonalization_inplace_selector<MatrixType,1,IsComplex>
507 {
508   typedef typename MatrixType::Scalar Scalar;
509 
510   template<typename DiagonalType, typename SubDiagonalType>
511   static void run(MatrixType& mat, DiagonalType& diag, SubDiagonalType&, bool extractQ)
512   {
513     diag(0,0) = numext::real(mat(0,0));
514     if(extractQ)
515       mat(0,0) = Scalar(1);
516   }
517 };
518 
519 /** \internal
520   * \eigenvalues_module \ingroup Eigenvalues_Module
521   *
522   * \brief Expression type for return value of Tridiagonalization::matrixT()
523   *
524   * \tparam MatrixType type of underlying dense matrix
525   */
526 template<typename MatrixType> struct TridiagonalizationMatrixTReturnType
527 : public ReturnByValue<TridiagonalizationMatrixTReturnType<MatrixType> >
528 {
529   public:
530     /** \brief Constructor.
531       *
532       * \param[in] mat The underlying dense matrix
533       */
534     TridiagonalizationMatrixTReturnType(const MatrixType& mat) : m_matrix(mat) { }
535 
536     template <typename ResultType>
537     inline void evalTo(ResultType& result) const
538     {
539       result.setZero();
540       result.template diagonal<1>() = m_matrix.template diagonal<-1>().conjugate();
541       result.diagonal() = m_matrix.diagonal();
542       result.template diagonal<-1>() = m_matrix.template diagonal<-1>();
543     }
544 
545     Index rows() const { return m_matrix.rows(); }
546     Index cols() const { return m_matrix.cols(); }
547 
548   protected:
549     typename MatrixType::Nested m_matrix;
550 };
551 
552 } // end namespace internal
553 
554 } // end namespace Eigen
555 
556 #endif // EIGEN_TRIDIAGONALIZATION_H
557