Eigen::LDLT< MatrixType_, UpLo_ > Class Template Reference

Robust Cholesky decomposition of a matrix with pivoting. More...

+ Inheritance diagram for Eigen::LDLT< MatrixType_, UpLo_ >:

Public Types

enum  {
  MaxRowsAtCompileTime ,
  MaxColsAtCompileTime ,
  UpLo
}
 
typedef SolverBase< LDLTBase
 
typedef MatrixType_ MatrixType
 
typedef PermutationMatrix< RowsAtCompileTime, MaxRowsAtCompileTimePermutationType
 
typedef Matrix< Scalar, RowsAtCompileTime, 1, 0, MaxRowsAtCompileTime, 1 > TmpMatrixType
 
typedef internal::LDLT_Traits< MatrixType, UpLoTraits
 
typedef Transpositions< RowsAtCompileTime, MaxRowsAtCompileTimeTranspositionType
 
- Public Types inherited from Eigen::SolverBase< LDLT< MatrixType_, UpLo_ > >
enum  
 
typedef std::conditional_t< NumTraits< Scalar >::IsComplex, CwiseUnaryOp< internal::scalar_conjugate_op< Scalar >, const ConstTransposeReturnType >, const ConstTransposeReturnTypeAdjointReturnType
 
typedef EigenBase< LDLT< MatrixType_, UpLo_ > > Base
 
typedef Scalar CoeffReturnType
 
typedef Transpose< const LDLT< MatrixType_, UpLo_ > > ConstTransposeReturnType
 
typedef internal::traits< LDLT< MatrixType_, UpLo_ > >::Scalar Scalar
 
- Public Types inherited from Eigen::EigenBase< Derived >
typedef Eigen::Index Index
 The interface type of indices. More...
 
typedef internal::traits< Derived >::StorageKind StorageKind
 

Public Member Functions

const LDLTadjoint () const
 
EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
template<typename InputType >
LDLT< MatrixType, UpLo_ > & compute (const EigenBase< InputType > &a)
 
template<typename InputType >
LDLTcompute (const EigenBase< InputType > &matrix)
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
bool isNegative (void) const
 
bool isPositive () const
 
 LDLT ()
 Default Constructor. More...
 
template<typename InputType >
 LDLT (const EigenBase< InputType > &matrix)
 Constructor with decomposition. More...
 
template<typename InputType >
 LDLT (EigenBase< InputType > &matrix)
 Constructs a LDLT factorization from a given matrix. More...
 
 LDLT (Index size)
 Default Constructor with memory preallocation. More...
 
Traits::MatrixL matrixL () const
 
const MatrixTypematrixLDLT () const
 
Traits::MatrixU matrixU () const
 
template<typename Derived >
LDLTrankUpdate (const MatrixBase< Derived > &w, const RealScalar &alpha=1)
 
template<typename Derived >
LDLT< MatrixType, UpLo_ > & rankUpdate (const MatrixBase< Derived > &w, const typename LDLT< MatrixType, UpLo_ >::RealScalar &sigma)
 
RealScalar rcond () const
 
MatrixType reconstructedMatrix () const
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
void setZero ()
 
template<typename Rhs >
const Solve< LDLT, Rhs > solve (const MatrixBase< Rhs > &b) const
 
template<typename Derived >
bool solveInPlace (MatrixBase< Derived > &bAndX) const
 
const TranspositionTypetranspositionsP () const
 
Diagonal< const MatrixTypevectorD () const
 
- Public Member Functions inherited from Eigen::SolverBase< LDLT< MatrixType_, UpLo_ > >
const AdjointReturnType adjoint () const
 
LDLT< MatrixType_, UpLo_ > & derived ()
 
const LDLT< MatrixType_, UpLo_ > & derived () const
 
const Solve< LDLT< MatrixType_, UpLo_ >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
 SolverBase ()
 
const ConstTransposeReturnType transpose () const
 
 ~SolverBase ()
 
- Public Member Functions inherited from Eigen::EigenBase< Derived >
template<typename Dest >
void addTo (Dest &dst) const
 
template<typename Dest >
void applyThisOnTheLeft (Dest &dst) const
 
template<typename Dest >
void applyThisOnTheRight (Dest &dst) const
 
EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
Derived & const_cast_derived () const
 
const Derived & const_derived () const
 
Derived & derived ()
 
const Derived & derived () const
 
template<typename Dest >
void evalTo (Dest &dst) const
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
EIGEN_CONSTEXPR Index size () const EIGEN_NOEXCEPT
 
template<typename Dest >
void subTo (Dest &dst) const
 

Protected Attributes

ComputationInfo m_info
 
bool m_isInitialized
 
RealScalar m_l1_norm
 
MatrixType m_matrix
 
internal::SignMatrix m_sign
 
TmpMatrixType m_temporary
 
TranspositionType m_transpositions
 

Additional Inherited Members

- Protected Member Functions inherited from Eigen::SolverBase< LDLT< MatrixType_, UpLo_ > >
void _check_solve_assertion (const Rhs &b) const
 

Detailed Description

template<typename MatrixType_, int UpLo_>
class Eigen::LDLT< MatrixType_, UpLo_ >

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 decomposition: Lower (default) or Upper. The other triangular part won't be read.

Perform a robust Cholesky decomposition of a positive semidefinite or negative semidefinite matrix \( A \) such that \( A = P^TLDL^*P \), 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 D 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.

See also
MatrixBase::ldlt(), SelfAdjointView::ldlt(), class LLT

Definition at line 61 of file LDLT.h.

Member Typedef Documentation

◆ Base

template<typename MatrixType_ , int UpLo_>
typedef SolverBase<LDLT> Eigen::LDLT< MatrixType_, UpLo_ >::Base

Definition at line 66 of file LDLT.h.

◆ MatrixType

template<typename MatrixType_ , int UpLo_>
typedef MatrixType_ Eigen::LDLT< MatrixType_, UpLo_ >::MatrixType

Definition at line 65 of file LDLT.h.

◆ PermutationType

template<typename MatrixType_ , int UpLo_>
typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime> Eigen::LDLT< MatrixType_, UpLo_ >::PermutationType

Definition at line 78 of file LDLT.h.

◆ TmpMatrixType

template<typename MatrixType_ , int UpLo_>
typedef Matrix<Scalar, RowsAtCompileTime, 1, 0, MaxRowsAtCompileTime, 1> Eigen::LDLT< MatrixType_, UpLo_ >::TmpMatrixType

Definition at line 75 of file LDLT.h.

◆ Traits

template<typename MatrixType_ , int UpLo_>
typedef internal::LDLT_Traits<MatrixType,UpLo> Eigen::LDLT< MatrixType_, UpLo_ >::Traits

Definition at line 80 of file LDLT.h.

◆ TranspositionType

template<typename MatrixType_ , int UpLo_>
typedef Transpositions<RowsAtCompileTime, MaxRowsAtCompileTime> Eigen::LDLT< MatrixType_, UpLo_ >::TranspositionType

Definition at line 77 of file LDLT.h.

Member Enumeration Documentation

◆ anonymous enum

template<typename MatrixType_ , int UpLo_>
anonymous enum
Enumerator
MaxRowsAtCompileTime 
MaxColsAtCompileTime 
UpLo 

Definition at line 70 of file LDLT.h.

70  {
71  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
72  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
73  UpLo = UpLo_
74  };
@ MaxColsAtCompileTime
Definition: LDLT.h:72
@ MaxRowsAtCompileTime
Definition: LDLT.h:71
@ UpLo
Definition: LDLT.h:73

Constructor & Destructor Documentation

◆ LDLT() [1/4]

template<typename MatrixType_ , int UpLo_>
Eigen::LDLT< MatrixType_, UpLo_ >::LDLT ( )
inline

Default Constructor.

The default constructor is useful in cases in which the user intends to perform decompositions via LDLT::compute(const MatrixType&).

Definition at line 87 of file LDLT.h.

88  : m_matrix(),
91  m_isInitialized(false)
92  {}
internal::SignMatrix m_sign
Definition: LDLT.h:287
bool m_isInitialized
Definition: LDLT.h:288
TranspositionType m_transpositions
Definition: LDLT.h:285
MatrixType m_matrix
Definition: LDLT.h:283

◆ LDLT() [2/4]

template<typename MatrixType_ , int UpLo_>
Eigen::LDLT< MatrixType_, UpLo_ >::LDLT ( Index  size)
inlineexplicit

Default Constructor with memory preallocation.

Like the default constructor but with preallocation of the internal data according to the specified problem size.

See also
LDLT()

Definition at line 100 of file LDLT.h.

101  : m_matrix(size, size),
103  m_temporary(size),
105  m_isInitialized(false)
106  {}
TmpMatrixType m_temporary
Definition: LDLT.h:286
EIGEN_CONSTEXPR Index size() const EIGEN_NOEXCEPT
Definition: EigenBase.h:69

◆ LDLT() [3/4]

template<typename MatrixType_ , int UpLo_>
template<typename InputType >
Eigen::LDLT< MatrixType_, UpLo_ >::LDLT ( const EigenBase< InputType > &  matrix)
inlineexplicit

Constructor with decomposition.

This calculates the decomposition for the input matrix.

See also
LDLT(Index size)

Definition at line 115 of file LDLT.h.

116  : m_matrix(matrix.rows(), matrix.cols()),
117  m_transpositions(matrix.rows()),
118  m_temporary(matrix.rows()),
120  m_isInitialized(false)
121  {
122  compute(matrix.derived());
123  }
LDLT & compute(const EigenBase< InputType > &matrix)

◆ LDLT() [4/4]

template<typename MatrixType_ , int UpLo_>
template<typename InputType >
Eigen::LDLT< MatrixType_, UpLo_ >::LDLT ( EigenBase< InputType > &  matrix)
inlineexplicit

Constructs a LDLT factorization from a given matrix.

This overloaded constructor is provided for inplace decomposition when MatrixType is a Eigen::Ref.

See also
LDLT(const EigenBase&)

Definition at line 132 of file LDLT.h.

133  : m_matrix(matrix.derived()),
134  m_transpositions(matrix.rows()),
135  m_temporary(matrix.rows()),
137  m_isInitialized(false)
138  {
139  compute(matrix.derived());
140  }

Member Function Documentation

◆ adjoint()

template<typename MatrixType_ , int UpLo_>
const LDLT& Eigen::LDLT< MatrixType_, UpLo_ >::adjoint ( ) const
inline
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)
Array< int, 3, 1 > b

Definition at line 249 of file LDLT.h.

249 { return *this; }

◆ cols()

template<typename MatrixType_ , int UpLo_>
EIGEN_CONSTEXPR Index Eigen::LDLT< MatrixType_, UpLo_ >::cols ( ) const
inline

Definition at line 252 of file LDLT.h.

252 { return m_matrix.cols(); }

◆ compute() [1/2]

template<typename MatrixType_ , int UpLo_>
template<typename InputType >
LDLT<MatrixType,UpLo_>& Eigen::LDLT< MatrixType_, UpLo_ >::compute ( const EigenBase< InputType > &  a)

Compute / recompute the LDLT decomposition A = L D L^* = U^* D U of matrix

Definition at line 498 of file LDLT.h.

499 {
500  eigen_assert(a.rows()==a.cols());
501  const Index size = a.rows();
502 
503  m_matrix = a.derived();
504 
505  // Compute matrix L1 norm = max abs column sum.
506  m_l1_norm = RealScalar(0);
507  // TODO move this code to SelfAdjointView
508  for (Index col = 0; col < size; ++col) {
509  RealScalar abs_col_sum;
510  if (UpLo_ == Lower)
511  abs_col_sum = m_matrix.col(col).tail(size - col).template lpNorm<1>() + m_matrix.row(col).head(col).template lpNorm<1>();
512  else
513  abs_col_sum = m_matrix.col(col).head(col).template lpNorm<1>() + m_matrix.row(col).tail(size - col).template lpNorm<1>();
514  if (abs_col_sum > m_l1_norm)
515  m_l1_norm = abs_col_sum;
516  }
517 
519  m_isInitialized = false;
522 
523  m_info = internal::ldlt_inplace<UpLo>::unblocked(m_matrix, m_transpositions, m_temporary, m_sign) ? Success : NumericalIssue;
524 
525  m_isInitialized = true;
526  return *this;
527 }
ColXpr col(Index i)
This is the const version of col().
#define eigen_assert(x)
Definition: Macros.h:902
ComputationInfo m_info
Definition: LDLT.h:289
RealScalar m_l1_norm
Definition: LDLT.h:284
constexpr void resize(Index rows, Index cols)
void resize(Index newSize)
@ Lower
Definition: Constants.h:211
@ NumericalIssue
Definition: Constants.h:448
@ Success
Definition: Constants.h:446
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:41

◆ compute() [2/2]

template<typename MatrixType_ , int UpLo_>
template<typename InputType >
LDLT& Eigen::LDLT< MatrixType_, UpLo_ >::compute ( const EigenBase< InputType > &  matrix)

◆ info()

template<typename MatrixType_ , int UpLo_>
ComputationInfo Eigen::LDLT< MatrixType_, UpLo_ >::info ( ) const
inline

Reports whether previous computation was successful.

Returns
Success if computation was successful, NumericalIssue if the factorization failed because of a zero pivot.

Definition at line 259 of file LDLT.h.

260  {
261  eigen_assert(m_isInitialized && "LDLT is not initialized.");
262  return m_info;
263  }

◆ isNegative()

template<typename MatrixType_ , int UpLo_>
bool Eigen::LDLT< MatrixType_, UpLo_ >::isNegative ( void  ) const
inline
Returns
true if the matrix is negative (semidefinite)

Definition at line 187 of file LDLT.h.

188  {
189  eigen_assert(m_isInitialized && "LDLT is not initialized.");
191  }
@ NegativeSemiDef
Definition: LDLT.h:33

◆ isPositive()

template<typename MatrixType_ , int UpLo_>
bool Eigen::LDLT< MatrixType_, UpLo_ >::isPositive ( ) const
inline
Returns
true if the matrix is positive (semidefinite)

Definition at line 180 of file LDLT.h.

181  {
182  eigen_assert(m_isInitialized && "LDLT is not initialized.");
184  }
@ PositiveSemiDef
Definition: LDLT.h:33

◆ matrixL()

template<typename MatrixType_ , int UpLo_>
Traits::MatrixL Eigen::LDLT< MatrixType_, UpLo_ >::matrixL ( ) const
inline
Returns
a view of the lower triangular matrix L

Definition at line 158 of file LDLT.h.

159  {
160  eigen_assert(m_isInitialized && "LDLT is not initialized.");
161  return Traits::getL(m_matrix);
162  }

◆ matrixLDLT()

template<typename MatrixType_ , int UpLo_>
const MatrixType& Eigen::LDLT< MatrixType_, UpLo_ >::matrixLDLT ( ) const
inline
Returns
the internal LDLT decomposition matrix

TODO: document the storage layout

Definition at line 236 of file LDLT.h.

237  {
238  eigen_assert(m_isInitialized && "LDLT is not initialized.");
239  return m_matrix;
240  }

◆ matrixU()

template<typename MatrixType_ , int UpLo_>
Traits::MatrixU Eigen::LDLT< MatrixType_, UpLo_ >::matrixU ( ) const
inline
Returns
a view of the upper triangular matrix U

Definition at line 151 of file LDLT.h.

152  {
153  eigen_assert(m_isInitialized && "LDLT is not initialized.");
154  return Traits::getU(m_matrix);
155  }

◆ rankUpdate() [1/2]

template<typename MatrixType_ , int UpLo_>
template<typename Derived >
LDLT& Eigen::LDLT< MatrixType_, UpLo_ >::rankUpdate ( const MatrixBase< Derived > &  w,
const RealScalar &  alpha = 1 
)

◆ rankUpdate() [2/2]

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 
)

Update the LDLT decomposition: given A = L D L^T, efficiently compute the decomposition of A + sigma w w^T.

Parameters
wa vector to be incorporated into the decomposition.
sigmaa scalar, +1 for updates and -1 for "downdates," which correspond to removing previously-added column vectors. Optional; default value is +1.
See also
setZero()

Definition at line 536 of file LDLT.h.

537 {
538  typedef typename TranspositionType::StorageIndex IndexType;
539  const Index size = w.rows();
540  if (m_isInitialized)
541  {
542  eigen_assert(m_matrix.rows()==size);
543  }
544  else
545  {
546  m_matrix.resize(size,size);
547  m_matrix.setZero();
549  for (Index i = 0; i < size; i++)
550  m_transpositions.coeffRef(i) = IndexType(i);
553  m_isInitialized = true;
554  }
555 
556  internal::ldlt_inplace<UpLo>::update(m_matrix, m_transpositions, m_temporary, w, sigma);
557 
558  return *this;
559 }
RowVector3d w
StorageIndex & coeffRef(Index i)

◆ rcond()

template<typename MatrixType_ , int UpLo_>
RealScalar Eigen::LDLT< MatrixType_, UpLo_ >::rcond ( ) const
inline
Returns
an estimate of the reciprocal condition number of the matrix of which *this is the LDLT decomposition.

Definition at line 223 of file LDLT.h.

224  {
225  eigen_assert(m_isInitialized && "LDLT is not initialized.");
227  }
Decomposition::RealScalar rcond_estimate_helper(typename Decomposition::RealScalar matrix_norm, const Decomposition &dec)
Reciprocal condition number estimator.

◆ reconstructedMatrix()

template<typename MatrixType , int UpLo_>
MatrixType Eigen::LDLT< MatrixType, UpLo_ >::reconstructedMatrix
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.

Definition at line 640 of file LDLT.h.

641 {
642  eigen_assert(m_isInitialized && "LDLT is not initialized.");
643  const Index size = m_matrix.rows();
645 
646  // P
647  res.setIdentity();
648  res = transpositionsP() * res;
649  // L^* P
650  res = matrixU() * res;
651  // D(L^*P)
652  res = vectorD().real().asDiagonal() * res;
653  // L(DL^*P)
654  res = matrixL() * res;
655  // P^T (LDL^*P)
657 
658  return res;
659 }
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
Matrix< float, 1, Dynamic > MatrixType
const TranspositionType & transpositionsP() const
Definition: LDLT.h:166
Traits::MatrixL matrixL() const
Definition: LDLT.h:158
Diagonal< const MatrixType > vectorD() const
Definition: LDLT.h:173
Traits::MatrixU matrixU() const
Definition: LDLT.h:151
Transpose< TranspositionsBase > transpose() const

◆ rows()

template<typename MatrixType_ , int UpLo_>
EIGEN_CONSTEXPR Index Eigen::LDLT< MatrixType_, UpLo_ >::rows ( ) const
inline

Definition at line 251 of file LDLT.h.

251 { return m_matrix.rows(); }

◆ setZero()

template<typename MatrixType_ , int UpLo_>
void Eigen::LDLT< MatrixType_, UpLo_ >::setZero ( )
inline

Clear any existing decomposition

See also
rankUpdate(w,sigma)

Definition at line 145 of file LDLT.h.

146  {
147  m_isInitialized = false;
148  }

◆ solve()

template<typename MatrixType_ , int UpLo_>
template<typename Rhs >
const Solve<LDLT, Rhs> Eigen::LDLT< MatrixType_, UpLo_ >::solve ( const MatrixBase< Rhs > &  b) const
inline
Returns
a solution x of \( A x = b \) 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::isApprox() directly, for instance like this:

bool a_solution_exists = (A*result).isApprox(b, precision);
MatrixXcf A
bool isApprox(const Scalar &x, const Scalar &y, const typename NumTraits< Scalar >::Real &precision=NumTraits< Scalar >::dummy_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 \( A x = b \) using the decomposition \( A = P^T L D L^* P \) by solving the systems \( P^T y_1 = b \), \( L y_2 = y_1 \), \( D y_3 = y_2 \), \( L^* y_4 = y_3 \) and \( P x = y_4 \) in succession. If the matrix \( A \) is singular, then \( D \) will also be singular (all the other matrices are invertible). In that case, the least-square solution of \( D y_3 = y_2 \) is computed. This does not mean that this function computes the least-square solution of \( A x = b \) if \( A \) is singular.

See also
MatrixBase::ldlt(), SelfAdjointView::ldlt()

◆ solveInPlace()

template<typename MatrixType , int UpLo_>
template<typename Derived >
bool Eigen::LDLT< MatrixType, UpLo_ >::solveInPlace ( MatrixBase< Derived > &  bAndX) const

Definition at line 626 of file LDLT.h.

627 {
628  eigen_assert(m_isInitialized && "LDLT is not initialized.");
629  eigen_assert(m_matrix.rows() == bAndX.rows());
630 
631  bAndX = this->solve(bAndX);
632 
633  return true;
634 }
const Solve< LDLT, Rhs > solve(const MatrixBase< Rhs > &b) const

◆ transpositionsP()

template<typename MatrixType_ , int UpLo_>
const TranspositionType& Eigen::LDLT< MatrixType_, UpLo_ >::transpositionsP ( ) const
inline
Returns
the permutation matrix P as a transposition sequence.

Definition at line 166 of file LDLT.h.

167  {
168  eigen_assert(m_isInitialized && "LDLT is not initialized.");
169  return m_transpositions;
170  }

◆ vectorD()

template<typename MatrixType_ , int UpLo_>
Diagonal<const MatrixType> Eigen::LDLT< MatrixType_, UpLo_ >::vectorD ( ) const
inline
Returns
the coefficients of the diagonal matrix D

Definition at line 173 of file LDLT.h.

174  {
175  eigen_assert(m_isInitialized && "LDLT is not initialized.");
176  return m_matrix.diagonal();
177  }

Member Data Documentation

◆ m_info

template<typename MatrixType_ , int UpLo_>
ComputationInfo Eigen::LDLT< MatrixType_, UpLo_ >::m_info
protected

Definition at line 289 of file LDLT.h.

◆ m_isInitialized

template<typename MatrixType_ , int UpLo_>
bool Eigen::LDLT< MatrixType_, UpLo_ >::m_isInitialized
protected

Definition at line 288 of file LDLT.h.

◆ m_l1_norm

template<typename MatrixType_ , int UpLo_>
RealScalar Eigen::LDLT< MatrixType_, UpLo_ >::m_l1_norm
protected

Definition at line 284 of file LDLT.h.

◆ m_matrix

template<typename MatrixType_ , int UpLo_>
MatrixType Eigen::LDLT< MatrixType_, UpLo_ >::m_matrix
protected

Definition at line 283 of file LDLT.h.

◆ m_sign

template<typename MatrixType_ , int UpLo_>
internal::SignMatrix Eigen::LDLT< MatrixType_, UpLo_ >::m_sign
protected

Definition at line 287 of file LDLT.h.

◆ m_temporary

template<typename MatrixType_ , int UpLo_>
TmpMatrixType Eigen::LDLT< MatrixType_, UpLo_ >::m_temporary
protected

Definition at line 286 of file LDLT.h.

◆ m_transpositions

template<typename MatrixType_ , int UpLo_>
TranspositionType Eigen::LDLT< MatrixType_, UpLo_ >::m_transpositions
protected

Definition at line 285 of file LDLT.h.


The documentation for this class was generated from the following file: