1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2008-2011 Gael Guennebaud <gael.guennebaud@inria.fr>
5// Copyright (C) 2009 Keir Mierle <mierle@gmail.com>
6// Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
7// Copyright (C) 2011 Timothy E. Holy <tim.holy@gmail.com >
8//
9// This Source Code Form is subject to the terms of the Mozilla
10// Public License v. 2.0. If a copy of the MPL was not distributed
11// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
12
13#ifndef EIGEN_LDLT_H
14#define EIGEN_LDLT_H
15
16namespace Eigen {
17
18namespace internal {
19  template<typename MatrixType, int UpLo> struct LDLT_Traits;
20
21  // PositiveSemiDef means positive semi-definite and non-zero; same for NegativeSemiDef
22  enum SignMatrix { PositiveSemiDef, NegativeSemiDef, ZeroSign, Indefinite };
23}
24
25/** \ingroup Cholesky_Module
26  *
27  * \class LDLT
28  *
29  * \brief Robust Cholesky decomposition of a matrix with pivoting
30  *
31  * \param MatrixType the type of the matrix of which to compute the LDL^T Cholesky decomposition
32  * \param UpLo the triangular part that will be used for the decompositon: Lower (default) or Upper.
33  *             The other triangular part won't be read.
34  *
35  * Perform a robust Cholesky decomposition of a positive semidefinite or negative semidefinite
36  * matrix \f$ A \f$ such that \f$ A =  P^TLDL^*P \f$, where P is a permutation matrix, L
37  * is lower triangular with a unit diagonal and D is a diagonal matrix.
38  *
39  * The decomposition uses pivoting to ensure stability, so that L will have
40  * zeros in the bottom right rank(A) - n submatrix. Avoiding the square root
41  * on D also stabilizes the computation.
42  *
43  * Remember that Cholesky decompositions are not rank-revealing. Also, do not use a Cholesky
44  * decomposition to determine whether a system of equations has a solution.
45  *
46  * \sa MatrixBase::ldlt(), class LLT
47  */
48template<typename _MatrixType, int _UpLo> class LDLT
49{
50  public:
51    typedef _MatrixType MatrixType;
52    enum {
53      RowsAtCompileTime = MatrixType::RowsAtCompileTime,
54      ColsAtCompileTime = MatrixType::ColsAtCompileTime,
55      Options = MatrixType::Options & ~RowMajorBit, // these are the options for the TmpMatrixType, we need a ColMajor matrix here!
56      MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
57      MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
58      UpLo = _UpLo
59    };
60    typedef typename MatrixType::Scalar Scalar;
61    typedef typename NumTraits<typename MatrixType::Scalar>::Real RealScalar;
62    typedef typename MatrixType::Index Index;
63    typedef Matrix<Scalar, RowsAtCompileTime, 1, Options, MaxRowsAtCompileTime, 1> TmpMatrixType;
64
65    typedef Transpositions<RowsAtCompileTime, MaxRowsAtCompileTime> TranspositionType;
66    typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime> PermutationType;
67
68    typedef internal::LDLT_Traits<MatrixType,UpLo> Traits;
69
70    /** \brief Default Constructor.
71      *
72      * The default constructor is useful in cases in which the user intends to
73      * perform decompositions via LDLT::compute(const MatrixType&).
74      */
75    LDLT()
76      : m_matrix(),
77        m_transpositions(),
78        m_sign(internal::ZeroSign),
79        m_isInitialized(false)
80    {}
81
82    /** \brief Default Constructor with memory preallocation
83      *
84      * Like the default constructor but with preallocation of the internal data
85      * according to the specified problem \a size.
86      * \sa LDLT()
87      */
88    LDLT(Index size)
89      : m_matrix(size, size),
90        m_transpositions(size),
91        m_temporary(size),
92        m_sign(internal::ZeroSign),
93        m_isInitialized(false)
94    {}
95
96    /** \brief Constructor with decomposition
97      *
98      * This calculates the decomposition for the input \a matrix.
99      * \sa LDLT(Index size)
100      */
101    LDLT(const MatrixType& matrix)
102      : m_matrix(matrix.rows(), matrix.cols()),
103        m_transpositions(matrix.rows()),
104        m_temporary(matrix.rows()),
105        m_sign(internal::ZeroSign),
106        m_isInitialized(false)
107    {
108      compute(matrix);
109    }
110
111    /** Clear any existing decomposition
112     * \sa rankUpdate(w,sigma)
113     */
114    void setZero()
115    {
116      m_isInitialized = false;
117    }
118
119    /** \returns a view of the upper triangular matrix U */
120    inline typename Traits::MatrixU matrixU() const
121    {
122      eigen_assert(m_isInitialized && "LDLT is not initialized.");
123      return Traits::getU(m_matrix);
124    }
125
126    /** \returns a view of the lower triangular matrix L */
127    inline typename Traits::MatrixL matrixL() const
128    {
129      eigen_assert(m_isInitialized && "LDLT is not initialized.");
130      return Traits::getL(m_matrix);
131    }
132
133    /** \returns the permutation matrix P as a transposition sequence.
134      */
135    inline const TranspositionType& transpositionsP() const
136    {
137      eigen_assert(m_isInitialized && "LDLT is not initialized.");
138      return m_transpositions;
139    }
140
141    /** \returns the coefficients of the diagonal matrix D */
142    inline Diagonal<const MatrixType> vectorD() const
143    {
144      eigen_assert(m_isInitialized && "LDLT is not initialized.");
145      return m_matrix.diagonal();
146    }
147
148    /** \returns true if the matrix is positive (semidefinite) */
149    inline bool isPositive() const
150    {
151      eigen_assert(m_isInitialized && "LDLT is not initialized.");
152      return m_sign == internal::PositiveSemiDef || m_sign == internal::ZeroSign;
153    }
154
155    #ifdef EIGEN2_SUPPORT
156    inline bool isPositiveDefinite() const
157    {
158      return isPositive();
159    }
160    #endif
161
162    /** \returns true if the matrix is negative (semidefinite) */
163    inline bool isNegative(void) const
164    {
165      eigen_assert(m_isInitialized && "LDLT is not initialized.");
166      return m_sign == internal::NegativeSemiDef || m_sign == internal::ZeroSign;
167    }
168
169    /** \returns a solution x of \f$ A x = b \f$ using the current decomposition of A.
170      *
171      * This function also supports in-place solves using the syntax <tt>x = decompositionObject.solve(x)</tt> .
172      *
173      * \note_about_checking_solutions
174      *
175      * More precisely, this method solves \f$ A x = b \f$ using the decomposition \f$ A = P^T L D L^* P \f$
176      * by solving the systems \f$ P^T y_1 = b \f$, \f$ L y_2 = y_1 \f$, \f$ D y_3 = y_2 \f$,
177      * \f$ L^* y_4 = y_3 \f$ and \f$ P x = y_4 \f$ in succession. If the matrix \f$ A \f$ is singular, then
178      * \f$ D \f$ will also be singular (all the other matrices are invertible). In that case, the
179      * least-square solution of \f$ D y_3 = y_2 \f$ is computed. This does not mean that this function
180      * computes the least-square solution of \f$ A x = b \f$ is \f$ A \f$ is singular.
181      *
182      * \sa MatrixBase::ldlt()
183      */
184    template<typename Rhs>
185    inline const internal::solve_retval<LDLT, Rhs>
186    solve(const MatrixBase<Rhs>& b) const
187    {
188      eigen_assert(m_isInitialized && "LDLT is not initialized.");
189      eigen_assert(m_matrix.rows()==b.rows()
190                && "LDLT::solve(): invalid number of rows of the right hand side matrix b");
191      return internal::solve_retval<LDLT, Rhs>(*this, b.derived());
192    }
193
194    #ifdef EIGEN2_SUPPORT
195    template<typename OtherDerived, typename ResultType>
196    bool solve(const MatrixBase<OtherDerived>& b, ResultType *result) const
197    {
198      *result = this->solve(b);
199      return true;
200    }
201    #endif
202
203    template<typename Derived>
204    bool solveInPlace(MatrixBase<Derived> &bAndX) const;
205
206    LDLT& compute(const MatrixType& matrix);
207
208    template <typename Derived>
209    LDLT& rankUpdate(const MatrixBase<Derived>& w, const RealScalar& alpha=1);
210
211    /** \returns the internal LDLT decomposition matrix
212      *
213      * TODO: document the storage layout
214      */
215    inline const MatrixType& matrixLDLT() const
216    {
217      eigen_assert(m_isInitialized && "LDLT is not initialized.");
218      return m_matrix;
219    }
220
221    MatrixType reconstructedMatrix() const;
222
223    inline Index rows() const { return m_matrix.rows(); }
224    inline Index cols() const { return m_matrix.cols(); }
225
226    /** \brief Reports whether previous computation was successful.
227      *
228      * \returns \c Success if computation was succesful,
229      *          \c NumericalIssue if the matrix.appears to be negative.
230      */
231    ComputationInfo info() const
232    {
233      eigen_assert(m_isInitialized && "LDLT is not initialized.");
234      return Success;
235    }
236
237  protected:
238
239    /** \internal
240      * Used to compute and store the Cholesky decomposition A = L D L^* = U^* D U.
241      * The strict upper part is used during the decomposition, the strict lower
242      * part correspond to the coefficients of L (its diagonal is equal to 1 and
243      * is not stored), and the diagonal entries correspond to D.
244      */
245    MatrixType m_matrix;
246    TranspositionType m_transpositions;
247    TmpMatrixType m_temporary;
248    internal::SignMatrix m_sign;
249    bool m_isInitialized;
250};
251
252namespace internal {
253
254template<int UpLo> struct ldlt_inplace;
255
256template<> struct ldlt_inplace<Lower>
257{
258  template<typename MatrixType, typename TranspositionType, typename Workspace>
259  static bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign)
260  {
261    using std::abs;
262    typedef typename MatrixType::Scalar Scalar;
263    typedef typename MatrixType::RealScalar RealScalar;
264    typedef typename MatrixType::Index Index;
265    eigen_assert(mat.rows()==mat.cols());
266    const Index size = mat.rows();
267
268    if (size <= 1)
269    {
270      transpositions.setIdentity();
271      if (numext::real(mat.coeff(0,0)) > 0) sign = PositiveSemiDef;
272      else if (numext::real(mat.coeff(0,0)) < 0) sign = NegativeSemiDef;
273      else sign = ZeroSign;
274      return true;
275    }
276
277    for (Index k = 0; k < size; ++k)
278    {
279      // Find largest diagonal element
280      Index index_of_biggest_in_corner;
281      mat.diagonal().tail(size-k).cwiseAbs().maxCoeff(&index_of_biggest_in_corner);
282      index_of_biggest_in_corner += k;
283
284      transpositions.coeffRef(k) = index_of_biggest_in_corner;
285      if(k != index_of_biggest_in_corner)
286      {
287        // apply the transposition while taking care to consider only
288        // the lower triangular part
289        Index s = size-index_of_biggest_in_corner-1; // trailing size after the biggest element
290        mat.row(k).head(k).swap(mat.row(index_of_biggest_in_corner).head(k));
291        mat.col(k).tail(s).swap(mat.col(index_of_biggest_in_corner).tail(s));
292        std::swap(mat.coeffRef(k,k),mat.coeffRef(index_of_biggest_in_corner,index_of_biggest_in_corner));
293        for(int i=k+1;i<index_of_biggest_in_corner;++i)
294        {
295          Scalar tmp = mat.coeffRef(i,k);
296          mat.coeffRef(i,k) = numext::conj(mat.coeffRef(index_of_biggest_in_corner,i));
297          mat.coeffRef(index_of_biggest_in_corner,i) = numext::conj(tmp);
298        }
299        if(NumTraits<Scalar>::IsComplex)
300          mat.coeffRef(index_of_biggest_in_corner,k) = numext::conj(mat.coeff(index_of_biggest_in_corner,k));
301      }
302
303      // partition the matrix:
304      //       A00 |  -  |  -
305      // lu  = A10 | A11 |  -
306      //       A20 | A21 | A22
307      Index rs = size - k - 1;
308      Block<MatrixType,Dynamic,1> A21(mat,k+1,k,rs,1);
309      Block<MatrixType,1,Dynamic> A10(mat,k,0,1,k);
310      Block<MatrixType,Dynamic,Dynamic> A20(mat,k+1,0,rs,k);
311
312      if(k>0)
313      {
314        temp.head(k) = mat.diagonal().real().head(k).asDiagonal() * A10.adjoint();
315        mat.coeffRef(k,k) -= (A10 * temp.head(k)).value();
316        if(rs>0)
317          A21.noalias() -= A20 * temp.head(k);
318      }
319
320      // In some previous versions of Eigen (e.g., 3.2.1), the scaling was omitted if the pivot
321      // was smaller than the cutoff value. However, soince LDLT is not rank-revealing
322      // we should only make sure we do not introduce INF or NaN values.
323      // LAPACK also uses 0 as the cutoff value.
324      RealScalar realAkk = numext::real(mat.coeffRef(k,k));
325      if((rs>0) && (abs(realAkk) > RealScalar(0)))
326        A21 /= realAkk;
327
328      if (sign == PositiveSemiDef) {
329        if (realAkk < 0) sign = Indefinite;
330      } else if (sign == NegativeSemiDef) {
331        if (realAkk > 0) sign = Indefinite;
332      } else if (sign == ZeroSign) {
333        if (realAkk > 0) sign = PositiveSemiDef;
334        else if (realAkk < 0) sign = NegativeSemiDef;
335      }
336    }
337
338    return true;
339  }
340
341  // Reference for the algorithm: Davis and Hager, "Multiple Rank
342  // Modifications of a Sparse Cholesky Factorization" (Algorithm 1)
343  // Trivial rearrangements of their computations (Timothy E. Holy)
344  // allow their algorithm to work for rank-1 updates even if the
345  // original matrix is not of full rank.
346  // Here only rank-1 updates are implemented, to reduce the
347  // requirement for intermediate storage and improve accuracy
348  template<typename MatrixType, typename WDerived>
349  static bool updateInPlace(MatrixType& mat, MatrixBase<WDerived>& w, const typename MatrixType::RealScalar& sigma=1)
350  {
351    using numext::isfinite;
352    typedef typename MatrixType::Scalar Scalar;
353    typedef typename MatrixType::RealScalar RealScalar;
354    typedef typename MatrixType::Index Index;
355
356    const Index size = mat.rows();
357    eigen_assert(mat.cols() == size && w.size()==size);
358
359    RealScalar alpha = 1;
360
361    // Apply the update
362    for (Index j = 0; j < size; j++)
363    {
364      // Check for termination due to an original decomposition of low-rank
365      if (!(isfinite)(alpha))
366        break;
367
368      // Update the diagonal terms
369      RealScalar dj = numext::real(mat.coeff(j,j));
370      Scalar wj = w.coeff(j);
371      RealScalar swj2 = sigma*numext::abs2(wj);
372      RealScalar gamma = dj*alpha + swj2;
373
374      mat.coeffRef(j,j) += swj2/alpha;
375      alpha += swj2/dj;
376
377
378      // Update the terms of L
379      Index rs = size-j-1;
380      w.tail(rs) -= wj * mat.col(j).tail(rs);
381      if(gamma != 0)
382        mat.col(j).tail(rs) += (sigma*numext::conj(wj)/gamma)*w.tail(rs);
383    }
384    return true;
385  }
386
387  template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
388  static bool update(MatrixType& mat, const TranspositionType& transpositions, Workspace& tmp, const WType& w, const typename MatrixType::RealScalar& sigma=1)
389  {
390    // Apply the permutation to the input w
391    tmp = transpositions * w;
392
393    return ldlt_inplace<Lower>::updateInPlace(mat,tmp,sigma);
394  }
395};
396
397template<> struct ldlt_inplace<Upper>
398{
399  template<typename MatrixType, typename TranspositionType, typename Workspace>
400  static EIGEN_STRONG_INLINE bool unblocked(MatrixType& mat, TranspositionType& transpositions, Workspace& temp, SignMatrix& sign)
401  {
402    Transpose<MatrixType> matt(mat);
403    return ldlt_inplace<Lower>::unblocked(matt, transpositions, temp, sign);
404  }
405
406  template<typename MatrixType, typename TranspositionType, typename Workspace, typename WType>
407  static EIGEN_STRONG_INLINE bool update(MatrixType& mat, TranspositionType& transpositions, Workspace& tmp, WType& w, const typename MatrixType::RealScalar& sigma=1)
408  {
409    Transpose<MatrixType> matt(mat);
410    return ldlt_inplace<Lower>::update(matt, transpositions, tmp, w.conjugate(), sigma);
411  }
412};
413
414template<typename MatrixType> struct LDLT_Traits<MatrixType,Lower>
415{
416  typedef const TriangularView<const MatrixType, UnitLower> MatrixL;
417  typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitUpper> MatrixU;
418  static inline MatrixL getL(const MatrixType& m) { return m; }
419  static inline MatrixU getU(const MatrixType& m) { return m.adjoint(); }
420};
421
422template<typename MatrixType> struct LDLT_Traits<MatrixType,Upper>
423{
424  typedef const TriangularView<const typename MatrixType::AdjointReturnType, UnitLower> MatrixL;
425  typedef const TriangularView<const MatrixType, UnitUpper> MatrixU;
426  static inline MatrixL getL(const MatrixType& m) { return m.adjoint(); }
427  static inline MatrixU getU(const MatrixType& m) { return m; }
428};
429
430} // end namespace internal
431
432/** Compute / recompute the LDLT decomposition A = L D L^* = U^* D U of \a matrix
433  */
434template<typename MatrixType, int _UpLo>
435LDLT<MatrixType,_UpLo>& LDLT<MatrixType,_UpLo>::compute(const MatrixType& a)
436{
437  eigen_assert(a.rows()==a.cols());
438  const Index size = a.rows();
439
440  m_matrix = a;
441
442  m_transpositions.resize(size);
443  m_isInitialized = false;
444  m_temporary.resize(size);
445
446  internal::ldlt_inplace<UpLo>::unblocked(m_matrix, m_transpositions, m_temporary, m_sign);
447
448  m_isInitialized = true;
449  return *this;
450}
451
452/** Update the LDLT decomposition:  given A = L D L^T, efficiently compute the decomposition of A + sigma w w^T.
453 * \param w a vector to be incorporated into the decomposition.
454 * \param sigma a scalar, +1 for updates and -1 for "downdates," which correspond to removing previously-added column vectors. Optional; default value is +1.
455 * \sa setZero()
456  */
457template<typename MatrixType, int _UpLo>
458template<typename Derived>
459LDLT<MatrixType,_UpLo>& LDLT<MatrixType,_UpLo>::rankUpdate(const MatrixBase<Derived>& w, const typename NumTraits<typename MatrixType::Scalar>::Real& sigma)
460{
461  const Index size = w.rows();
462  if (m_isInitialized)
463  {
464    eigen_assert(m_matrix.rows()==size);
465  }
466  else
467  {
468    m_matrix.resize(size,size);
469    m_matrix.setZero();
470    m_transpositions.resize(size);
471    for (Index i = 0; i < size; i++)
472      m_transpositions.coeffRef(i) = i;
473    m_temporary.resize(size);
474    m_sign = sigma>=0 ? internal::PositiveSemiDef : internal::NegativeSemiDef;
475    m_isInitialized = true;
476  }
477
478  internal::ldlt_inplace<UpLo>::update(m_matrix, m_transpositions, m_temporary, w, sigma);
479
480  return *this;
481}
482
483namespace internal {
484template<typename _MatrixType, int _UpLo, typename Rhs>
485struct solve_retval<LDLT<_MatrixType,_UpLo>, Rhs>
486  : solve_retval_base<LDLT<_MatrixType,_UpLo>, Rhs>
487{
488  typedef LDLT<_MatrixType,_UpLo> LDLTType;
489  EIGEN_MAKE_SOLVE_HELPERS(LDLTType,Rhs)
490
491  template<typename Dest> void evalTo(Dest& dst) const
492  {
493    eigen_assert(rhs().rows() == dec().matrixLDLT().rows());
494    // dst = P b
495    dst = dec().transpositionsP() * rhs();
496
497    // dst = L^-1 (P b)
498    dec().matrixL().solveInPlace(dst);
499
500    // dst = D^-1 (L^-1 P b)
501    // more precisely, use pseudo-inverse of D (see bug 241)
502    using std::abs;
503    using std::max;
504    typedef typename LDLTType::MatrixType MatrixType;
505    typedef typename LDLTType::Scalar Scalar;
506    typedef typename LDLTType::RealScalar RealScalar;
507    const typename Diagonal<const MatrixType>::RealReturnType vectorD(dec().vectorD());
508    // In some previous versions, tolerance was set to the max of 1/highest and the maximal diagonal entry * epsilon
509    // as motivated by LAPACK's xGELSS:
510    // RealScalar tolerance = (max)(vectorD.array().abs().maxCoeff() *NumTraits<RealScalar>::epsilon(),RealScalar(1) / NumTraits<RealScalar>::highest());
511    // However, LDLT is not rank revealing, and so adjusting the tolerance wrt to the highest
512    // diagonal element is not well justified and to numerical issues in some cases.
513    // Moreover, Lapack's xSYTRS routines use 0 for the tolerance.
514    RealScalar tolerance = RealScalar(1) / NumTraits<RealScalar>::highest();
515
516    for (Index i = 0; i < vectorD.size(); ++i) {
517      if(abs(vectorD(i)) > tolerance)
518        dst.row(i) /= vectorD(i);
519      else
520        dst.row(i).setZero();
521    }
522
523    // dst = L^-T (D^-1 L^-1 P b)
524    dec().matrixU().solveInPlace(dst);
525
526    // dst = P^-1 (L^-T D^-1 L^-1 P b) = A^-1 b
527    dst = dec().transpositionsP().transpose() * dst;
528  }
529};
530}
531
532/** \internal use x = ldlt_object.solve(x);
533  *
534  * This is the \em in-place version of solve().
535  *
536  * \param bAndX represents both the right-hand side matrix b and result x.
537  *
538  * \returns true always! If you need to check for existence of solutions, use another decomposition like LU, QR, or SVD.
539  *
540  * This version avoids a copy when the right hand side matrix b is not
541  * needed anymore.
542  *
543  * \sa LDLT::solve(), MatrixBase::ldlt()
544  */
545template<typename MatrixType,int _UpLo>
546template<typename Derived>
547bool LDLT<MatrixType,_UpLo>::solveInPlace(MatrixBase<Derived> &bAndX) const
548{
549  eigen_assert(m_isInitialized && "LDLT is not initialized.");
550  eigen_assert(m_matrix.rows() == bAndX.rows());
551
552  bAndX = this->solve(bAndX);
553
554  return true;
555}
556
557/** \returns the matrix represented by the decomposition,
558 * i.e., it returns the product: P^T L D L^* P.
559 * This function is provided for debug purpose. */
560template<typename MatrixType, int _UpLo>
561MatrixType LDLT<MatrixType,_UpLo>::reconstructedMatrix() const
562{
563  eigen_assert(m_isInitialized && "LDLT is not initialized.");
564  const Index size = m_matrix.rows();
565  MatrixType res(size,size);
566
567  // P
568  res.setIdentity();
569  res = transpositionsP() * res;
570  // L^* P
571  res = matrixU() * res;
572  // D(L^*P)
573  res = vectorD().real().asDiagonal() * res;
574  // L(DL^*P)
575  res = matrixL() * res;
576  // P^T (LDL^*P)
577  res = transpositionsP().transpose() * res;
578
579  return res;
580}
581
582/** \cholesky_module
583  * \returns the Cholesky decomposition with full pivoting without square root of \c *this
584  */
585template<typename MatrixType, unsigned int UpLo>
586inline const LDLT<typename SelfAdjointView<MatrixType, UpLo>::PlainObject, UpLo>
587SelfAdjointView<MatrixType, UpLo>::ldlt() const
588{
589  return LDLT<PlainObject,UpLo>(m_matrix);
590}
591
592/** \cholesky_module
593  * \returns the Cholesky decomposition with full pivoting without square root of \c *this
594  */
595template<typename Derived>
596inline const LDLT<typename MatrixBase<Derived>::PlainObject>
597MatrixBase<Derived>::ldlt() const
598{
599  return LDLT<PlainObject>(derived());
600}
601
602} // end namespace Eigen
603
604#endif // EIGEN_LDLT_H
605