Standard Cholesky decomposition (LL^T) of a matrix and associated features. More...
Public Types | |
enum | { MaxColsAtCompileTime } |
enum | { PacketSize , AlignmentMask , UpLo } |
typedef SolverBase< LLT > | Base |
typedef MatrixType_ | MatrixType |
typedef internal::LLT_Traits< MatrixType, UpLo > | Traits |
Public Types inherited from Eigen::SolverBase< LLT< MatrixType_, UpLo_ > > | |
enum | |
typedef std::conditional_t< NumTraits< Scalar >::IsComplex, CwiseUnaryOp< internal::scalar_conjugate_op< Scalar >, const ConstTransposeReturnType >, const ConstTransposeReturnType > | AdjointReturnType |
typedef EigenBase< LLT< MatrixType_, UpLo_ > > | Base |
typedef Scalar | CoeffReturnType |
typedef Transpose< const LLT< MatrixType_, UpLo_ > > | ConstTransposeReturnType |
typedef internal::traits< LLT< 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 LLT & | adjoint () const EIGEN_NOEXCEPT |
EIGEN_CONSTEXPR Index | cols () const EIGEN_NOEXCEPT |
template<typename InputType > | |
LLT< MatrixType, UpLo_ > & | compute (const EigenBase< InputType > &a) |
template<typename InputType > | |
LLT & | compute (const EigenBase< InputType > &matrix) |
ComputationInfo | info () const |
Reports whether previous computation was successful. More... | |
LLT () | |
Default Constructor. More... | |
template<typename InputType > | |
LLT (const EigenBase< InputType > &matrix) | |
template<typename InputType > | |
LLT (EigenBase< InputType > &matrix) | |
Constructs a LLT factorization from a given matrix. More... | |
LLT (Index size) | |
Default Constructor with memory preallocation. More... | |
Traits::MatrixL | matrixL () const |
const MatrixType & | matrixLLT () const |
Traits::MatrixU | matrixU () const |
template<typename VectorType > | |
LLT< MatrixType_, UpLo_ > & | rankUpdate (const VectorType &v, const RealScalar &sigma) |
template<typename VectorType > | |
LLT & | rankUpdate (const VectorType &vec, const RealScalar &sigma=1) |
RealScalar | rcond () const |
MatrixType | reconstructedMatrix () const |
EIGEN_CONSTEXPR Index | rows () const EIGEN_NOEXCEPT |
template<typename Rhs > | |
const Solve< LLT, Rhs > | solve (const MatrixBase< Rhs > &b) const |
template<typename Derived > | |
void | solveInPlace (const MatrixBase< Derived > &bAndX) const |
Public Member Functions inherited from Eigen::SolverBase< LLT< MatrixType_, UpLo_ > > | |
const AdjointReturnType | adjoint () const |
LLT< MatrixType_, UpLo_ > & | derived () |
const LLT< MatrixType_, UpLo_ > & | derived () const |
const Solve< LLT< 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 |
Additional Inherited Members | |
Protected Member Functions inherited from Eigen::SolverBase< LLT< MatrixType_, UpLo_ > > | |
void | _check_solve_assertion (const Rhs &b) const |
Standard Cholesky decomposition (LL^T) of a matrix and associated features.
MatrixType_ | the type of the matrix of which we are computing the LL^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. |
This class performs a LL^T Cholesky decomposition of a symmetric, positive definite matrix A such that A = LL^* = U^*U, where L is lower triangular.
While the Cholesky decomposition is particularly useful to solve selfadjoint problems like D^*D x = b, for that purpose, we recommend the Cholesky decomposition without square root which is more stable and even faster. Nevertheless, this standard Cholesky decomposition remains useful in many other situations like generalised eigen problems with hermitian matrices.
Remember that Cholesky decompositions are not rank-revealing. This LLT decomposition is only stable on positive definite matrices, use LDLT instead for the semidefinite case. Also, do not use a Cholesky decomposition to determine whether a system of equations has a solution.
Example:
Output:
The matrix A is 4 -1 2 -1 6 0 2 0 5 The Cholesky factor L is 2 0 0 -0.5 2.4 0 1 0.209 1.99 To check this, let us compute L * L.transpose() 4 -1 2 -1 6 0 2 0 5 This should equal the matrix A
Performance: for best performance, it is recommended to use a column-major storage format with the Lower triangular part (the default), or, equivalently, a row-major storage format with the Upper triangular part. Otherwise, you might get a 20% slowdown for the full factorization step, and rank-updates can be up to 3 times slower.
This class supports the inplace decomposition mechanism.
Note that during the decomposition, only the lower (or upper, as defined by UpLo_) triangular part of A is considered. Therefore, the strict lower part does not have to store correct values.
typedef SolverBase<LLT> Eigen::LLT< MatrixType_, UpLo_ >::Base |
typedef MatrixType_ Eigen::LLT< MatrixType_, UpLo_ >::MatrixType |
typedef internal::LLT_Traits<MatrixType,UpLo> Eigen::LLT< MatrixType_, UpLo_ >::Traits |
anonymous enum |
anonymous enum |
|
inline |
|
inlineexplicit |
|
inlineexplicit |
|
inlineexplicit |
Constructs a LLT factorization from a given matrix.
This overloaded constructor is provided for inplace decomposition when MatrixType
is a Eigen::Ref.
Definition at line 122 of file LLT.h.
|
inline |
*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:
|
inline |
LLT<MatrixType,UpLo_>& Eigen::LLT< MatrixType_, UpLo_ >::compute | ( | const EigenBase< InputType > & | a | ) |
Computes / recomputes the Cholesky decomposition A = LL^* = U^*U of matrix
Example:
Output:
Here is the matrix A: 2 -1 -1 3 Here is the right hand side b: 1 2 3 1 Computing LLT decomposition... The solution is: 1.2 1.4 1.4 0.8 The matrix A is now: 2 -1 -1 4 Computing LLT decomposition... The solution is now: 1 1.29 1 0.571
Definition at line 431 of file LLT.h.
LLT& Eigen::LLT< MatrixType_, UpLo_ >::compute | ( | const EigenBase< InputType > & | matrix | ) |
|
inline |
Reports whether previous computation was successful.
Success
if computation was successful, NumericalIssue
if the matrix.appears not to be positive definite.
|
inline |
|
inline |
TODO: document the storage layout
|
inline |
LLT<MatrixType_,UpLo_>& Eigen::LLT< MatrixType_, UpLo_ >::rankUpdate | ( | const VectorType & | v, |
const RealScalar & | sigma | ||
) |
Performs a rank one update (or dowdate) of the current decomposition. If A = LL^* before the rank one update, then after it we have LL^* = A + sigma * v v^* where v must be a vector of same dimension.
Definition at line 466 of file LLT.h.
LLT& Eigen::LLT< MatrixType_, UpLo_ >::rankUpdate | ( | const VectorType & | vec, |
const RealScalar & | sigma = 1 |
||
) |
|
inline |
*this
is the Cholesky decomposition. Definition at line 168 of file LLT.h.
MatrixType Eigen::LLT< MatrixType, UpLo_ >::reconstructedMatrix |
|
inline |
|
inline |
Since this LLT class assumes anyway that the matrix A is invertible, the solution theoretically exists and is unique regardless of b.
Example:
Output:
2.02 2.97
void Eigen::LLT< MatrixType, UpLo_ >::solveInPlace | ( | const MatrixBase< Derived > & | bAndX | ) | const |
|
protected |
|
protected |
|
protected |
|
protected |