template<typename _MatrixType, int _UpLo>
          LDLT class
        
        Robust Cholesky decomposition of a matrix with pivoting.
| Template parameters | |
|---|---|
| _MatrixType | the type of the matrix of which to compute the LDL^T Cholesky decomposition | 
| _UpLo | the triangular part that will be used for the decompositon: Lower (default) or Upper. The other triangular part won't be read. | 
Contents
Perform a robust Cholesky decomposition of a positive semidefinite or negative semidefinite matrix such that , where P is a permutation matrix, L is lower triangular with a unit diagonal and D is a diagonal matrix.
The decomposition uses pivoting to ensure stability, so that L will have zeros in the bottom right rank(A) - n submatrix. Avoiding the square root on D also stabilizes the computation.
Remember that Cholesky decompositions are not rank-revealing. Also, do not use a Cholesky decomposition to determine whether a system of equations has a solution.
This class supports the inplace decomposition mechanism.
Base classes
- 
              template<typename Derived>class SolverBase
 - A base class for matrix decomposition and solvers.
 
Constructors, destructors, conversion operators
- LDLT()
 - Default Constructor.
 - LDLT(Index size) explicit
 - Default Constructor with memory preallocation.
 - 
              template<typename InputType>LDLT(const EigenBase<InputType>& matrix) explicit
 - Constructor with decomposition.
 - 
              template<typename InputType>LDLT(EigenBase<InputType>& matrix) explicit
 - Constructs a LDLT factorization from a given matrix.
 
Public functions
- auto adjoint() const -> const LDLT&
 - 
              template<typename InputType>auto compute(const EigenBase<InputType>& a) -> LDLT<MatrixType, _UpLo>&
 - auto info() const -> ComputationInfo
 - Reports whether previous computation was successful.
 - auto isNegative(void) const -> bool
 - auto isPositive() const -> bool
 - auto matrixL() const -> Traits::MatrixL
 - auto matrixLDLT() const -> const MatrixType&
 - auto matrixU() const -> Traits::MatrixU
 - 
              template<typename Derived>auto rankUpdate(const MatrixBase<Derived>& w, const typename LDLT<MatrixType, _UpLo>::RealScalar& sigma) -> LDLT<MatrixType, _UpLo>&
 - auto rcond() const -> RealScalar
 - auto reconstructedMatrix() const -> MatrixType
 - void setZero()
 - 
              template<typename Rhs>auto solve(const MatrixBase<Rhs>& b) const -> const Solve<LDLT, Rhs>
 - auto transpositionsP() const -> const TranspositionType&
 - auto vectorD() const -> Diagonal<const MatrixType>
 
Function documentation
              
                template<typename _MatrixType, int _UpLo>
              
               Eigen:: LDLT<_MatrixType, _UpLo>:: LDLT()
            
            Default Constructor.
The default constructor is useful in cases in which the user intends to perform decompositions via LDLT::compute(const MatrixType&).
              
                template<typename _MatrixType, int _UpLo>
                template<typename InputType>
              
               Eigen:: LDLT<_MatrixType, _UpLo>:: LDLT(EigenBase<InputType>& matrix) explicit 
            
            Constructs a LDLT factorization from a given matrix.
This overloaded constructor is provided for inplace decomposition when MatrixType is a Eigen::
              
                template<typename _MatrixType, int _UpLo>
              
              const LDLT& Eigen:: LDLT<_MatrixType, _UpLo>:: adjoint() const
            
            | Returns | the adjoint of *this, that is, a const reference to the decomposition itself as the underlying matrix is self-adjoint. | 
                
|---|
This method is provided for compatibility with other matrix decompositions, thus enabling generic code such as: x = decomposition.adjoint().solve(b)
              
                template<typename _MatrixType, int _UpLo>
              
              ComputationInfo Eigen:: LDLT<_MatrixType, _UpLo>:: info() const
            
            Reports whether previous computation was successful.
| Returns | Success if computation was successful, NumericalIssue if the factorization failed because of a zero pivot. | 
                
|---|
              
                template<typename _MatrixType, int _UpLo>
              
              bool Eigen:: LDLT<_MatrixType, _UpLo>:: isNegative(void) const
            
            | Returns | true if the matrix is negative (semidefinite) | 
|---|
              
                template<typename _MatrixType, int _UpLo>
              
              bool Eigen:: LDLT<_MatrixType, _UpLo>:: isPositive() const
            
            | Returns | true if the matrix is positive (semidefinite) | 
|---|
              
                template<typename _MatrixType, int _UpLo>
              
              Traits::MatrixL Eigen:: LDLT<_MatrixType, _UpLo>:: matrixL() const
            
            | Returns | a view of the lower triangular matrix L | 
|---|
              
                template<typename _MatrixType, int _UpLo>
              
              const MatrixType& Eigen:: LDLT<_MatrixType, _UpLo>:: matrixLDLT() const
            
            | Returns | the internal LDLT decomposition matrix | 
|---|
TODO: document the storage layout
              
                template<typename _MatrixType, int _UpLo>
              
              Traits::MatrixU Eigen:: LDLT<_MatrixType, _UpLo>:: matrixU() const
            
            | Returns | a view of the upper triangular matrix U | 
|---|
              
                template<typename _MatrixType, int _UpLo>
                template<typename Derived>
              
              LDLT<MatrixType, _UpLo>& Eigen:: LDLT<_MatrixType, _UpLo>:: rankUpdate(const MatrixBase<Derived>& w,
              const typename LDLT<MatrixType, _UpLo>::RealScalar& sigma)
            
            | Parameters | |
|---|---|
| w | a vector to be incorporated into the decomposition. | 
| sigma | a scalar, +1 for updates and -1 for "downdates," which correspond to removing previously-added column vectors. Optional; default value is +1. | 
Update the LDLT decomposition: given A = L D L^T, efficiently compute the decomposition of A + sigma w w^T.
              
                template<typename _MatrixType, int _UpLo>
              
              MatrixType Eigen:: LDLT<_MatrixType, _UpLo>:: reconstructedMatrix() const
            
            | Returns | the matrix represented by the decomposition, i.e., it returns the product: P^T L D L^* P. This function is provided for debug purpose. | 
|---|
              
                template<typename _MatrixType, int _UpLo>
              
              void Eigen:: LDLT<_MatrixType, _UpLo>:: setZero()
            
Clear any existing decomposition
              
                template<typename _MatrixType, int _UpLo>
                template<typename Rhs>
              
              const Solve<LDLT, Rhs> Eigen:: LDLT<_MatrixType, _UpLo>:: solve(const MatrixBase<Rhs>& b) const
            
            | Returns | a solution x of using the current decomposition of A. | 
|---|
This function also supports in-place solves using the syntax x = decompositionObject.solve(x) .
This method just tries to find as good a solution as possible. If you want to check whether a solution exists or if it is accurate, just call this function to get a result and then compute the error of this result, or use MatrixBase::bool a_solution_exists = (A*result).isApprox(b, precision); This method avoids dividing by zero, so that the non-existence of a solution doesn't by itself mean that you'll get inf or nan values.
More precisely, this method solves using the decomposition by solving the systems , , , and in succession. If the matrix is singular, then will also be singular (all the other matrices are invertible). In that case, the least-square solution of is computed. This does not mean that this function computes the least-square solution of is is singular.
              
                template<typename _MatrixType, int _UpLo>
              
              const TranspositionType& Eigen:: LDLT<_MatrixType, _UpLo>:: transpositionsP() const
            
            | Returns | the permutation matrix P as a transposition sequence. | 
|---|