Prusa Slicer 2.6.0
Loading...
Searching...
No Matches
Eigen::LDLT< _MatrixType, _UpLo > Class Template Reference

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

#include <src/eigen/Eigen/src/Cholesky/LDLT.h>

+ Collaboration diagram for Eigen::LDLT< _MatrixType, _UpLo >:

Public Types

enum  {
  RowsAtCompileTime = MatrixType::RowsAtCompileTime , ColsAtCompileTime = MatrixType::ColsAtCompileTime , MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime , MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime ,
  UpLo = _UpLo
}
 
typedef _MatrixType MatrixType
 
typedef MatrixType::Scalar Scalar
 
typedef NumTraits< typenameMatrixType::Scalar >::Real RealScalar
 
typedef Eigen::Index Index
 
typedef MatrixType::StorageIndex StorageIndex
 
typedef Matrix< Scalar, RowsAtCompileTime, 1, 0, MaxRowsAtCompileTime, 1 > TmpMatrixType
 
typedef Transpositions< RowsAtCompileTime, MaxRowsAtCompileTimeTranspositionType
 
typedef PermutationMatrix< RowsAtCompileTime, MaxRowsAtCompileTimePermutationType
 
typedef internal::LDLT_Traits< MatrixType, UpLoTraits
 

Public Member Functions

 LDLT ()
 Default Constructor.
 
 LDLT (Index size)
 Default Constructor with memory preallocation.
 
template<typename InputType >
 LDLT (const EigenBase< InputType > &matrix)
 Constructor with decomposition.
 
template<typename InputType >
 LDLT (EigenBase< InputType > &matrix)
 Constructs a LDLT factorization from a given matrix.
 
void setZero ()
 
Traits::MatrixU matrixU () const
 
Traits::MatrixL matrixL () const
 
const TranspositionTypetranspositionsP () const
 
Diagonal< const MatrixTypevectorD () const
 
bool isPositive () const
 
bool isNegative (void) const
 
template<typename Rhs >
const Solve< LDLT, Rhs > solve (const MatrixBase< Rhs > &b) const
 
template<typename Derived >
bool solveInPlace (MatrixBase< Derived > &bAndX) const
 
template<typename InputType >
LDLTcompute (const EigenBase< InputType > &matrix)
 
RealScalar rcond () const
 
template<typename Derived >
LDLTrankUpdate (const MatrixBase< Derived > &w, const RealScalar &alpha=1)
 
const MatrixTypematrixLDLT () const
 
MatrixType reconstructedMatrix () const
 
const LDLTadjoint () const
 
Index rows () const
 
Index cols () const
 
ComputationInfo info () const
 Reports whether previous computation was successful.
 
template<typename RhsType , typename DstType >
EIGEN_DEVICE_FUNC void _solve_impl (const RhsType &rhs, DstType &dst) const
 
template<typename InputType >
LDLT< MatrixType, _UpLo > & compute (const EigenBase< InputType > &a)
 
template<typename Derived >
LDLT< MatrixType, _UpLo > & rankUpdate (const MatrixBase< Derived > &w, const typename LDLT< MatrixType, _UpLo >::RealScalar &sigma)
 
template<typename RhsType , typename DstType >
void _solve_impl (const RhsType &rhs, DstType &dst) const
 

Static Protected Member Functions

static void check_template_parameters ()
 

Protected Attributes

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

Detailed Description

template<typename _MatrixType, int _UpLo>
class Eigen::LDLT< _MatrixType, _UpLo >

Robust Cholesky decomposition of a matrix with pivoting.

Template Parameters
_MatrixTypethe type of the matrix of which to compute the LDL^T Cholesky decomposition
_UpLothe triangular part that will be used for the decompositon: 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 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.

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

Member Typedef Documentation

◆ Index

template<typename _MatrixType , int _UpLo>
typedef Eigen::Index Eigen::LDLT< _MatrixType, _UpLo >::Index

◆ MatrixType

template<typename _MatrixType , int _UpLo>
typedef _MatrixType Eigen::LDLT< _MatrixType, _UpLo >::MatrixType

◆ PermutationType

template<typename _MatrixType , int _UpLo>
typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime> Eigen::LDLT< _MatrixType, _UpLo >::PermutationType

◆ RealScalar

template<typename _MatrixType , int _UpLo>
typedef NumTraits<typenameMatrixType::Scalar>::Real Eigen::LDLT< _MatrixType, _UpLo >::RealScalar

◆ Scalar

template<typename _MatrixType , int _UpLo>
typedef MatrixType::Scalar Eigen::LDLT< _MatrixType, _UpLo >::Scalar

◆ StorageIndex

template<typename _MatrixType , int _UpLo>
typedef MatrixType::StorageIndex Eigen::LDLT< _MatrixType, _UpLo >::StorageIndex

◆ TmpMatrixType

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

◆ Traits

template<typename _MatrixType , int _UpLo>
typedef internal::LDLT_Traits<MatrixType,UpLo> Eigen::LDLT< _MatrixType, _UpLo >::Traits

◆ TranspositionType

template<typename _MatrixType , int _UpLo>
typedef Transpositions<RowsAtCompileTime, MaxRowsAtCompileTime> Eigen::LDLT< _MatrixType, _UpLo >::TranspositionType

Member Enumeration Documentation

◆ anonymous enum

template<typename _MatrixType , int _UpLo>
anonymous enum
Enumerator
RowsAtCompileTime 
ColsAtCompileTime 
MaxRowsAtCompileTime 
MaxColsAtCompileTime 
UpLo 
54 {
55 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
56 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
57 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
58 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime,
59 UpLo = _UpLo
60 };
@ MaxColsAtCompileTime
Definition LDLT.h:58
@ MaxRowsAtCompileTime
Definition LDLT.h:57
@ RowsAtCompileTime
Definition LDLT.h:55
@ UpLo
Definition LDLT.h:59
@ ColsAtCompileTime
Definition LDLT.h:56

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&).

78 : m_matrix(),
81 m_isInitialized(false)
82 {}
TranspositionType m_transpositions
Definition LDLT.h:280
MatrixType m_matrix
Definition LDLT.h:278
internal::SignMatrix m_sign
Definition LDLT.h:282
bool m_isInitialized
Definition LDLT.h:283
@ ZeroSign
Definition LDLT.h:22

◆ 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()
91 : m_matrix(size, size),
92 m_transpositions(size),
93 m_temporary(size),
95 m_isInitialized(false)
96 {}
TmpMatrixType m_temporary
Definition LDLT.h:281

◆ 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)
106 : m_matrix(matrix.rows(), matrix.cols()),
107 m_transpositions(matrix.rows()),
108 m_temporary(matrix.rows()),
110 m_isInitialized(false)
111 {
112 compute(matrix.derived());
113 }
LDLT & compute(const EigenBase< InputType > &matrix)

References Eigen::LDLT< _MatrixType, _UpLo >::compute(), and Eigen::EigenBase< Derived >::derived().

+ Here is the call graph for this function:

◆ 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&)
123 : m_matrix(matrix.derived()),
124 m_transpositions(matrix.rows()),
125 m_temporary(matrix.rows()),
127 m_isInitialized(false)
128 {
129 compute(matrix.derived());
130 }

References Eigen::LDLT< _MatrixType, _UpLo >::compute(), and Eigen::EigenBase< Derived >::derived().

+ Here is the call graph for this function:

Member Function Documentation

◆ _solve_impl() [1/2]

template<typename _MatrixType , int _UpLo>
template<typename RhsType , typename DstType >
EIGEN_DEVICE_FUNC void Eigen::LDLT< _MatrixType, _UpLo >::_solve_impl ( const RhsType &  rhs,
DstType &  dst 
) const

◆ _solve_impl() [2/2]

template<typename _MatrixType , int _UpLo>
template<typename RhsType , typename DstType >
void Eigen::LDLT< _MatrixType, _UpLo >::_solve_impl ( const RhsType &  rhs,
DstType &  dst 
) const
562{
563 eigen_assert(rhs.rows() == rows());
564 // dst = P b
565 dst = m_transpositions * rhs;
566
567 // dst = L^-1 (P b)
568 matrixL().solveInPlace(dst);
569
570 // dst = D^-1 (L^-1 P b)
571 // more precisely, use pseudo-inverse of D (see bug 241)
572 using std::abs;
573 const typename Diagonal<const MatrixType>::RealReturnType vecD(vectorD());
574 // In some previous versions, tolerance was set to the max of 1/highest (or rather numeric_limits::min())
575 // and the maximal diagonal entry * epsilon as motivated by LAPACK's xGELSS:
576 // RealScalar tolerance = numext::maxi(vecD.array().abs().maxCoeff() * NumTraits<RealScalar>::epsilon(),RealScalar(1) / NumTraits<RealScalar>::highest());
577 // However, LDLT is not rank revealing, and so adjusting the tolerance wrt to the highest
578 // diagonal element is not well justified and leads to numerical issues in some cases.
579 // Moreover, Lapack's xSYTRS routines use 0 for the tolerance.
580 // Using numeric_limits::min() gives us more robustness to denormals.
581 RealScalar tolerance = (std::numeric_limits<RealScalar>::min)();
582
583 for (Index i = 0; i < vecD.size(); ++i)
584 {
585 if(abs(vecD(i)) > tolerance)
586 dst.row(i) /= vecD(i);
587 else
588 dst.row(i).setZero();
589 }
590
591 // dst = L^-T (D^-1 L^-1 P b)
592 matrixU().solveInPlace(dst);
593
594 // dst = P^-1 (L^-T D^-1 L^-1 P b) = A^-1 b
595 dst = m_transpositions.transpose() * dst;
596}
#define eigen_assert(x)
Definition Macros.h:579
NumTraits< typenameMatrixType::Scalar >::Real RealScalar
Definition LDLT.h:62
Traits::MatrixU matrixU() const
Definition LDLT.h:141
Index rows() const
Definition LDLT.h:245
Diagonal< const MatrixType > vectorD() const
Definition LDLT.h:163
Eigen::Index Index
Definition LDLT.h:63
Traits::MatrixL matrixL() const
Definition LDLT.h:148
Transpose< TranspositionsBase > transpose() const
Definition Transpositions.h:112
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC half abs(const half &a)
Definition Half.h:445

References eigen_assert.

◆ 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)
243{ return *this; };

◆ check_template_parameters()

template<typename _MatrixType , int _UpLo>
static void Eigen::LDLT< _MatrixType, _UpLo >::check_template_parameters ( )
inlinestaticprotected
268 {
270 }
#define EIGEN_STATIC_ASSERT_NON_INTEGER(TYPE)
Definition StaticAssert.h:184
MatrixType::Scalar Scalar
Definition LDLT.h:61

References EIGEN_STATIC_ASSERT_NON_INTEGER.

◆ cols()

template<typename _MatrixType , int _UpLo>
Index Eigen::LDLT< _MatrixType, _UpLo >::cols ( ) const
inline
246{ return m_matrix.cols(); }

References Eigen::LDLT< _MatrixType, _UpLo >::m_matrix.

◆ 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

494{
496
497 eigen_assert(a.rows()==a.cols());
498 const Index size = a.rows();
499
500 m_matrix = a.derived();
501
502 // Compute matrix L1 norm = max abs column sum.
504 // TODO move this code to SelfAdjointView
505 for (Index col = 0; col < size; ++col) {
506 RealScalar abs_col_sum;
507 if (_UpLo == Lower)
508 abs_col_sum = m_matrix.col(col).tail(size - col).template lpNorm<1>() + m_matrix.row(col).head(col).template lpNorm<1>();
509 else
510 abs_col_sum = m_matrix.col(col).head(col).template lpNorm<1>() + m_matrix.row(col).tail(size - col).template lpNorm<1>();
511 if (abs_col_sum > m_l1_norm)
512 m_l1_norm = abs_col_sum;
513 }
514
516 m_isInitialized = false;
517 m_temporary.resize(size);
519
520 m_info = internal::ldlt_inplace<UpLo>::unblocked(m_matrix, m_transpositions, m_temporary, m_sign) ? Success : NumericalIssue;
521
522 m_isInitialized = true;
523 return *this;
524}
EIGEN_DEVICE_FUNC ColXpr col(Index i)
This is the const version of col().
Definition BlockMethods.h:838
ComputationInfo m_info
Definition LDLT.h:284
static void check_template_parameters()
Definition LDLT.h:267
RealScalar m_l1_norm
Definition LDLT.h:279
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
Definition PlainObjectBase.h:279
void resize(Index newSize)
Definition Transpositions.h:74
@ Lower
Definition Constants.h:204
@ NumericalIssue
Definition Constants.h:434
@ Success
Definition Constants.h:432
constexpr auto size(const C &c) -> decltype(c.size())
Definition span.hpp:183

References col(), eigen_assert, Eigen::Lower, Eigen::NumericalIssue, Eigen::Success, and Eigen::internal::ZeroSign.

+ Here is the call graph for this function:

◆ compute() [2/2]

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

Referenced by Eigen::LDLT< _MatrixType, _UpLo >::LDLT(), and Eigen::LDLT< _MatrixType, _UpLo >::LDLT().

+ Here is the caller graph for this function:

◆ 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 succesful, NumericalIssue if the factorization failed because of a zero pivot.
254 {
255 eigen_assert(m_isInitialized && "LDLT is not initialized.");
256 return m_info;
257 }

References eigen_assert, Eigen::LDLT< _MatrixType, _UpLo >::m_info, and Eigen::LDLT< _MatrixType, _UpLo >::m_isInitialized.

◆ isNegative()

template<typename _MatrixType , int _UpLo>
bool Eigen::LDLT< _MatrixType, _UpLo >::isNegative ( void  ) const
inline
Returns
true if the matrix is negative (semidefinite)
178 {
179 eigen_assert(m_isInitialized && "LDLT is not initialized.");
181 }
@ NegativeSemiDef
Definition LDLT.h:22

References eigen_assert, Eigen::LDLT< _MatrixType, _UpLo >::m_isInitialized, Eigen::LDLT< _MatrixType, _UpLo >::m_sign, Eigen::internal::NegativeSemiDef, and Eigen::internal::ZeroSign.

◆ isPositive()

template<typename _MatrixType , int _UpLo>
bool Eigen::LDLT< _MatrixType, _UpLo >::isPositive ( ) const
inline
Returns
true if the matrix is positive (semidefinite)
171 {
172 eigen_assert(m_isInitialized && "LDLT is not initialized.");
174 }
@ PositiveSemiDef
Definition LDLT.h:22

References eigen_assert, Eigen::LDLT< _MatrixType, _UpLo >::m_isInitialized, Eigen::LDLT< _MatrixType, _UpLo >::m_sign, Eigen::internal::PositiveSemiDef, and Eigen::internal::ZeroSign.

◆ matrixL()

template<typename _MatrixType , int _UpLo>
Traits::MatrixL Eigen::LDLT< _MatrixType, _UpLo >::matrixL ( ) const
inline
Returns
a view of the lower triangular matrix L
149 {
150 eigen_assert(m_isInitialized && "LDLT is not initialized.");
151 return Traits::getL(m_matrix);
152 }

References eigen_assert, Eigen::LDLT< _MatrixType, _UpLo >::m_isInitialized, and Eigen::LDLT< _MatrixType, _UpLo >::m_matrix.

◆ 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

231 {
232 eigen_assert(m_isInitialized && "LDLT is not initialized.");
233 return m_matrix;
234 }

References eigen_assert, Eigen::LDLT< _MatrixType, _UpLo >::m_isInitialized, and Eigen::LDLT< _MatrixType, _UpLo >::m_matrix.

◆ matrixU()

template<typename _MatrixType , int _UpLo>
Traits::MatrixU Eigen::LDLT< _MatrixType, _UpLo >::matrixU ( ) const
inline
Returns
a view of the upper triangular matrix U
142 {
143 eigen_assert(m_isInitialized && "LDLT is not initialized.");
144 return Traits::getU(m_matrix);
145 }

References eigen_assert, Eigen::LDLT< _MatrixType, _UpLo >::m_isInitialized, and Eigen::LDLT< _MatrixType, _UpLo >::m_matrix.

◆ 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()
534{
535 typedef typename TranspositionType::StorageIndex IndexType;
536 const Index size = w.rows();
537 if (m_isInitialized)
538 {
539 eigen_assert(m_matrix.rows()==size);
540 }
541 else
542 {
543 m_matrix.resize(size,size);
544 m_matrix.setZero();
546 for (Index i = 0; i < size; i++)
547 m_transpositions.coeffRef(i) = IndexType(i);
548 m_temporary.resize(size);
550 m_isInitialized = true;
551 }
552
553 internal::ldlt_inplace<UpLo>::update(m_matrix, m_transpositions, m_temporary, w, sigma);
554
555 return *this;
556}
StorageIndex & coeffRef(Index i)
Definition Transpositions.h:58
IndicesType::Scalar StorageIndex
Definition Transpositions.h:165

References eigen_assert, Eigen::internal::NegativeSemiDef, and Eigen::internal::PositiveSemiDef.

◆ 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.
218 {
219 eigen_assert(m_isInitialized && "LDLT is not initialized.");
221 }
Decomposition::RealScalar rcond_estimate_helper(typename Decomposition::RealScalar matrix_norm, const Decomposition &dec)
Reciprocal condition number estimator.
Definition ConditionEstimator.h:159

References eigen_assert, Eigen::LDLT< _MatrixType, _UpLo >::m_isInitialized, Eigen::LDLT< _MatrixType, _UpLo >::m_l1_norm, and Eigen::internal::rcond_estimate_helper().

+ Here is the call graph for this function:

◆ 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.
629{
630 eigen_assert(m_isInitialized && "LDLT is not initialized.");
631 const Index size = m_matrix.rows();
632 MatrixType res(size,size);
633
634 // P
635 res.setIdentity();
636 res = transpositionsP() * res;
637 // L^* P
638 res = matrixU() * res;
639 // D(L^*P)
640 res = vectorD().real().asDiagonal() * res;
641 // L(DL^*P)
642 res = matrixL() * res;
643 // P^T (LDL^*P)
644 res = transpositionsP().transpose() * res;
645
646 return res;
647}
_MatrixType MatrixType
Definition LDLT.h:53
const TranspositionType & transpositionsP() const
Definition LDLT.h:156

References eigen_assert.

◆ rows()

template<typename _MatrixType , int _UpLo>
Index Eigen::LDLT< _MatrixType, _UpLo >::rows ( ) const
inline
245{ return m_matrix.rows(); }

References Eigen::LDLT< _MatrixType, _UpLo >::m_matrix.

◆ setZero()

template<typename _MatrixType , int _UpLo>
void Eigen::LDLT< _MatrixType, _UpLo >::setZero ( )
inline

Clear any existing decomposition

See also
rankUpdate(w,sigma)
136 {
137 m_isInitialized = false;
138 }

References Eigen::LDLT< _MatrixType, _UpLo >::m_isInitialized.

◆ 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) .

\note_about_checking_solutions

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 $ is $ A $ is singular.

See also
MatrixBase::ldlt(), SelfAdjointView::ldlt()
201 {
202 eigen_assert(m_isInitialized && "LDLT is not initialized.");
203 eigen_assert(m_matrix.rows()==b.rows()
204 && "LDLT::solve(): invalid number of rows of the right hand side matrix b");
205 return Solve<LDLT, Rhs>(*this, b.derived());
206 }

References eigen_assert, Eigen::LDLT< _MatrixType, _UpLo >::m_isInitialized, and Eigen::LDLT< _MatrixType, _UpLo >::m_matrix.

◆ solveInPlace()

template<typename MatrixType , int _UpLo>
template<typename Derived >
bool Eigen::LDLT< MatrixType, _UpLo >::solveInPlace ( MatrixBase< Derived > &  bAndX) const
615{
616 eigen_assert(m_isInitialized && "LDLT is not initialized.");
617 eigen_assert(m_matrix.rows() == bAndX.rows());
618
619 bAndX = this->solve(bAndX);
620
621 return true;
622}
const Solve< LDLT, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition LDLT.h:200

References eigen_assert.

◆ transpositionsP()

template<typename _MatrixType , int _UpLo>
const TranspositionType & Eigen::LDLT< _MatrixType, _UpLo >::transpositionsP ( ) const
inline
Returns
the permutation matrix P as a transposition sequence.
157 {
158 eigen_assert(m_isInitialized && "LDLT is not initialized.");
159 return m_transpositions;
160 }

References eigen_assert, Eigen::LDLT< _MatrixType, _UpLo >::m_isInitialized, and Eigen::LDLT< _MatrixType, _UpLo >::m_transpositions.

◆ vectorD()

template<typename _MatrixType , int _UpLo>
Diagonal< const MatrixType > Eigen::LDLT< _MatrixType, _UpLo >::vectorD ( ) const
inline
Returns
the coefficients of the diagonal matrix D
164 {
165 eigen_assert(m_isInitialized && "LDLT is not initialized.");
166 return m_matrix.diagonal();
167 }

References eigen_assert, Eigen::LDLT< _MatrixType, _UpLo >::m_isInitialized, and Eigen::LDLT< _MatrixType, _UpLo >::m_matrix.

Member Data Documentation

◆ m_info

template<typename _MatrixType , int _UpLo>
ComputationInfo Eigen::LDLT< _MatrixType, _UpLo >::m_info
protected

◆ m_isInitialized

◆ m_l1_norm

template<typename _MatrixType , int _UpLo>
RealScalar Eigen::LDLT< _MatrixType, _UpLo >::m_l1_norm
protected

◆ m_matrix

◆ m_sign

template<typename _MatrixType , int _UpLo>
internal::SignMatrix Eigen::LDLT< _MatrixType, _UpLo >::m_sign
protected

◆ m_temporary

template<typename _MatrixType , int _UpLo>
TmpMatrixType Eigen::LDLT< _MatrixType, _UpLo >::m_temporary
protected

◆ m_transpositions

template<typename _MatrixType , int _UpLo>
TranspositionType Eigen::LDLT< _MatrixType, _UpLo >::m_transpositions
protected

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