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

Standard Cholesky decomposition (LL^T) of a matrix and associated features. More...

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

Public Types

enum  { MaxColsAtCompileTime }
 
enum  {
  PacketSize ,
  AlignmentMask ,
  UpLo
}
 
typedef SolverBase< LLTBase
 
typedef MatrixType_ MatrixType
 
typedef internal::LLT_Traits< MatrixType, UpLoTraits
 
- 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 ConstTransposeReturnTypeAdjointReturnType
 
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 LLTadjoint () const EIGEN_NOEXCEPT
 
EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
template<typename InputType >
LLT< MatrixType, UpLo_ > & compute (const EigenBase< InputType > &a)
 
template<typename InputType >
LLTcompute (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 MatrixTypematrixLLT () const
 
Traits::MatrixU matrixU () const
 
template<typename VectorType >
LLT< MatrixType_, UpLo_ > & rankUpdate (const VectorType &v, const RealScalar &sigma)
 
template<typename VectorType >
LLTrankUpdate (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
 

Detailed Description

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

Standard Cholesky decomposition (LL^T) of a matrix and associated features.

Template Parameters
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:

MatrixXd A(3,3);
A << 4,-1,2, -1,6,0, 2,0,5;
cout << "The matrix A is" << endl << A << endl;
LLT<MatrixXd> lltOfA(A); // compute the Cholesky decomposition of A
MatrixXd L = lltOfA.matrixL(); // retrieve factor L in the decomposition
// The previous two lines can also be written as "L = A.llt().matrixL()"
cout << "The Cholesky factor L is" << endl << L << endl;
cout << "To check this, let us compute L * L.transpose()" << endl;
cout << L * L.transpose() << endl;
cout << "This should equal the matrix A" << endl;
MatrixXcf A
MatrixXd L
Definition: LLT_example.cpp:6
A<< 4,-1, 2, -1, 6, 0, 2, 0, 5;cout<< "The matrix A is"<< endl<< A<< endl;LLT< MatrixXd > lltOfA(A)
Matrix< double, Dynamic, Dynamic > MatrixXd
DynamicĂ—Dynamic matrix of type double.
Definition: Matrix.h:502

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.

See also
MatrixBase::llt(), SelfAdjointView::llt(), class LDLT

Definition at line 68 of file LLT.h.

Member Typedef Documentation

◆ Base

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

Definition at line 73 of file LLT.h.

◆ MatrixType

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

Definition at line 72 of file LLT.h.

◆ Traits

template<typename MatrixType_ , int UpLo_>
typedef internal::LLT_Traits<MatrixType,UpLo> Eigen::LLT< MatrixType_, UpLo_ >::Traits

Definition at line 87 of file LLT.h.

Member Enumeration Documentation

◆ anonymous enum

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

Definition at line 77 of file LLT.h.

77  {
78  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
79  };
@ MaxColsAtCompileTime
Definition: LLT.h:78

◆ anonymous enum

template<typename MatrixType_ , int UpLo_>
anonymous enum
Enumerator
PacketSize 
AlignmentMask 
UpLo 

Definition at line 81 of file LLT.h.

81  {
83  AlignmentMask = int(PacketSize)-1,
84  UpLo = UpLo_
85  };
@ UpLo
Definition: LLT.h:84
@ PacketSize
Definition: LLT.h:82
@ AlignmentMask
Definition: LLT.h:83

Constructor & Destructor Documentation

◆ LLT() [1/4]

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

Default Constructor.

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

Definition at line 95 of file LLT.h.

95 : m_matrix(), m_isInitialized(false) {}
MatrixType m_matrix
Definition: LLT.h:228
bool m_isInitialized
Definition: LLT.h:230

◆ LLT() [2/4]

template<typename MatrixType_ , int UpLo_>
Eigen::LLT< MatrixType_, UpLo_ >::LLT ( 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
LLT()

Definition at line 103 of file LLT.h.

103  : m_matrix(size, size),
104  m_isInitialized(false) {}
EIGEN_CONSTEXPR Index size() const EIGEN_NOEXCEPT
Definition: EigenBase.h:69

◆ LLT() [3/4]

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

Definition at line 107 of file LLT.h.

108  : m_matrix(matrix.rows(), matrix.cols()),
109  m_isInitialized(false)
110  {
111  compute(matrix.derived());
112  }
LLT & compute(const EigenBase< InputType > &matrix)

◆ LLT() [4/4]

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

Constructs a LLT factorization from a given matrix.

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

See also
LLT(const EigenBase&)

Definition at line 122 of file LLT.h.

123  : m_matrix(matrix.derived()),
124  m_isInitialized(false)
125  {
126  compute(matrix.derived());
127  }

Member Function Documentation

◆ adjoint()

template<typename MatrixType_ , int UpLo_>
const LLT& Eigen::LLT< 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 204 of file LLT.h.

204 { return *this; }

◆ cols()

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

Definition at line 207 of file LLT.h.

207 { return m_matrix.cols(); }

◆ compute() [1/2]

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

Computes / recomputes the Cholesky decomposition A = LL^* = U^*U of matrix

Returns
a reference to *this

Example:

#include <iostream>
#include <Eigen/Dense>
int main()
{
A << 2, -1, -1, 3;
b << 1, 2, 3, 1;
std::cout << "Here is the matrix A:\n" << A << std::endl;
std::cout << "Here is the right hand side b:\n" << b << std::endl;
std::cout << "Computing LLT decomposition..." << std::endl;
llt.compute(A);
std::cout << "The solution is:\n" << llt.solve(b) << std::endl;
A(1,1)++;
std::cout << "The matrix A is now:\n" << A << std::endl;
std::cout << "Computing LLT decomposition..." << std::endl;
llt.compute(A);
std::cout << "The solution is now:\n" << llt.solve(b) << std::endl;
}
Standard Cholesky decomposition (LL^T) of a matrix and associated features.
Definition: LLT.h:70
const Solve< LLT, Rhs > solve(const MatrixBase< Rhs > &b) const
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:182
int main(int, char **)
Definition: class_Block.cpp:18

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.

432 {
433  eigen_assert(a.rows()==a.cols());
434  const Index size = a.rows();
435  m_matrix.resize(size, size);
436  if (!internal::is_same_dense(m_matrix, a.derived()))
437  m_matrix = a.derived();
438 
439  // Compute matrix L1 norm = max abs column sum.
440  m_l1_norm = RealScalar(0);
441  // TODO move this code to SelfAdjointView
442  for (Index col = 0; col < size; ++col) {
443  RealScalar abs_col_sum;
444  if (UpLo_ == Lower)
445  abs_col_sum = m_matrix.col(col).tail(size - col).template lpNorm<1>() + m_matrix.row(col).head(col).template lpNorm<1>();
446  else
447  abs_col_sum = m_matrix.col(col).head(col).template lpNorm<1>() + m_matrix.row(col).tail(size - col).template lpNorm<1>();
448  if (abs_col_sum > m_l1_norm)
449  m_l1_norm = abs_col_sum;
450  }
451 
452  m_isInitialized = true;
453  bool ok = Traits::inplace_decomposition(m_matrix);
454  m_info = ok ? Success : NumericalIssue;
455 
456  return *this;
457 }
ColXpr col(Index i)
This is the const version of col().
#define eigen_assert(x)
Definition: Macros.h:902
RealScalar m_l1_norm
Definition: LLT.h:229
ComputationInfo m_info
Definition: LLT.h:231
@ Lower
Definition: Constants.h:211
@ NumericalIssue
Definition: Constants.h:448
@ Success
Definition: Constants.h:446
bool is_same_dense(const T1 &mat1, const T2 &mat2, std::enable_if_t< possibly_same_dense< T1, T2 >::value > *=0)
Definition: XprHelper.h:745
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:41

◆ compute() [2/2]

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

◆ info()

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

Reports whether previous computation was successful.

Returns
Success if computation was successful, NumericalIssue if the matrix.appears not to be positive definite.

Definition at line 193 of file LLT.h.

194  {
195  eigen_assert(m_isInitialized && "LLT is not initialized.");
196  return m_info;
197  }

◆ matrixL()

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

Definition at line 137 of file LLT.h.

138  {
139  eigen_assert(m_isInitialized && "LLT is not initialized.");
140  return Traits::getL(m_matrix);
141  }

◆ matrixLLT()

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

TODO: document the storage layout

Definition at line 179 of file LLT.h.

180  {
181  eigen_assert(m_isInitialized && "LLT is not initialized.");
182  return m_matrix;
183  }

◆ matrixU()

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

Definition at line 130 of file LLT.h.

131  {
132  eigen_assert(m_isInitialized && "LLT is not initialized.");
133  return Traits::getU(m_matrix);
134  }

◆ rankUpdate() [1/2]

template<typename MatrixType_ , int UpLo_>
template<typename VectorType >
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.

467 {
469  eigen_assert(v.size()==m_matrix.cols());
471  if(internal::llt_inplace<typename MatrixType::Scalar, UpLo>::rankUpdate(m_matrix,v,sigma)>=0)
473  else
474  m_info = Success;
475 
476  return *this;
477 }
Array< int, Dynamic, 1 > v
#define EIGEN_STATIC_ASSERT_VECTOR_ONLY(TYPE)
Definition: StaticAssert.h:36

◆ rankUpdate() [2/2]

template<typename MatrixType_ , int UpLo_>
template<typename VectorType >
LLT& Eigen::LLT< MatrixType_, UpLo_ >::rankUpdate ( const VectorType &  vec,
const RealScalar &  sigma = 1 
)

◆ rcond()

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

Definition at line 168 of file LLT.h.

169  {
170  eigen_assert(m_isInitialized && "LLT is not initialized.");
171  eigen_assert(m_info == Success && "LLT failed because matrix appears to be negative");
173  }
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::LLT< MatrixType, UpLo_ >::reconstructedMatrix
Returns
the matrix represented by the decomposition, i.e., it returns the product: L L^*. This function is provided for debug purpose.

Definition at line 525 of file LLT.h.

526 {
527  eigen_assert(m_isInitialized && "LLT is not initialized.");
528  return matrixL() * matrixL().adjoint().toDenseMatrix();
529 }
Traits::MatrixL matrixL() const
Definition: LLT.h:137

◆ rows()

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

Definition at line 206 of file LLT.h.

206 { return m_matrix.rows(); }

◆ solve()

template<typename MatrixType_ , int UpLo_>
template<typename Rhs >
const Solve<LLT, Rhs> Eigen::LLT< MatrixType_, UpLo_ >::solve ( const MatrixBase< Rhs > &  b) const
inline
Returns
the solution x of \( A x = b \) using the current decomposition of A.

Since this LLT class assumes anyway that the matrix A is invertible, the solution theoretically exists and is unique regardless of b.

Example:

typedef Matrix<float,Dynamic,2> DataMatrix;
// let's generate some samples on the 3D plane of equation z = 2x+3y (with some noise)
DataMatrix samples = DataMatrix::Random(12,2);
VectorXf elevations = 2*samples.col(0) + 3*samples.col(1) + VectorXf::Random(12)*0.1;
// and let's solve samples * [x y]^T = elevations in least square sense:
Matrix<float,2,1> xy
= (samples.adjoint() * samples).llt().solve((samples.adjoint()*elevations));
cout << xy << endl;
Matrix< float, 2, 1 > xy
Definition: LLT_solve.cpp:7
Matrix< float, Dynamic, 2 > DataMatrix
Definition: LLT_solve.cpp:1
DataMatrix samples
Definition: LLT_solve.cpp:3
VectorXf elevations
Definition: LLT_solve.cpp:4
static const RandomReturnType Random()
Definition: Random.h:114
Matrix< float, Dynamic, 1 > VectorXf
DynamicĂ—1 vector of type float.
Definition: Matrix.h:501

Output:

2.02
2.97
See also
solveInPlace(), MatrixBase::llt(), SelfAdjointView::llt()

◆ solveInPlace()

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

Definition at line 513 of file LLT.h.

514 {
515  eigen_assert(m_isInitialized && "LLT is not initialized.");
516  eigen_assert(m_matrix.rows()==bAndX.rows());
517  matrixL().solveInPlace(bAndX);
518  matrixU().solveInPlace(bAndX);
519 }
Traits::MatrixU matrixU() const
Definition: LLT.h:130

Member Data Documentation

◆ m_info

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

Definition at line 231 of file LLT.h.

◆ m_isInitialized

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

Definition at line 230 of file LLT.h.

◆ m_l1_norm

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

Definition at line 229 of file LLT.h.

◆ m_matrix

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

Definition at line 228 of file LLT.h.


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