Eigen  3.2.7
 All Classes Namespaces Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
LLT< _MatrixType, _UpLo > Class Template Reference

Detailed Description

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

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

Parameters
MatrixTypethe type of the matrix of which we are computing the LL^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.

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;

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
See Also
MatrixBase::llt(), class LDLT

Public Member Functions

LLTcompute (const MatrixType &matrix)
 
ComputationInfo info () const
 Reports whether previous computation was successful. More...
 
 LLT ()
 Default Constructor. 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)
 
MatrixType reconstructedMatrix () const
 
template<typename Rhs >
const internal::solve_retval
< LLT, Rhs > 
solve (const MatrixBase< Rhs > &b) const
 

Constructor & Destructor Documentation

LLT ( )
inline

Default Constructor.

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

LLT ( Index  size)
inline

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()

Member Function Documentation

LLT< MatrixType, _UpLo > & compute ( const MatrixType &  a)

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

Returns
a reference to *this

Example:

#include <iostream>
#include <Eigen/Dense>
using namespace std;
using namespace Eigen;
int main()
{
Matrix2f A, b;
A << 2, -1, -1, 3;
b << 1, 2, 3, 1;
cout << "Here is the matrix A:\n" << A << endl;
cout << "Here is the right hand side b:\n" << b << endl;
cout << "Computing LLT decomposition..." << endl;
llt.compute(A);
cout << "The solution is:\n" << llt.solve(b) << endl;
A(1,1)++;
cout << "The matrix A is now:\n" << A << endl;
cout << "Computing LLT decomposition..." << endl;
llt.compute(A);
cout << "The solution is now:\n" << llt.solve(b) << endl;
}

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

References Eigen::NumericalIssue, and Eigen::Success.

ComputationInfo info ( ) const
inline

Reports whether previous computation was successful.

Returns
Success if computation was succesful, NumericalIssue if the matrix.appears to be negative.
Traits::MatrixL matrixL ( ) const
inline
Returns
a view of the lower triangular matrix L

Referenced by GeneralizedSelfAdjointEigenSolver< _MatrixType >::compute().

const MatrixType& matrixLLT ( ) const
inline
Returns
the LLT decomposition matrix

TODO: document the storage layout

Traits::MatrixU matrixU ( ) const
inline
Returns
a view of the upper triangular matrix U

Referenced by GeneralizedSelfAdjointEigenSolver< _MatrixType >::compute().

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.

References Eigen::NumericalIssue, and Eigen::Success.

MatrixType reconstructedMatrix ( ) const
Returns
the matrix represented by the decomposition, i.e., it returns the product: L L^*. This function is provided for debug purpose.
const internal::solve_retval<LLT, Rhs> 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;

Output:

2.02
2.97
See Also
solveInPlace(), MatrixBase::llt()

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