LU< MatrixType > Class Template Reference

LU decomposition of a matrix with complete pivoting, and related features. More...

#include <LU.h>

List of all members.

Public Types

enum  { MaxSmallDimAtCompileTime }
typedef Matrix< Scalar,
MatrixType::RowsAtCompileTime,
1, MatrixType::Options,
MatrixType::MaxRowsAtCompileTime, 1 > 
ColVectorType
typedef Matrix< typename
MatrixType::Scalar,
MatrixType::RowsAtCompileTime,
Dynamic, MatrixType::Options,
MatrixType::MaxRowsAtCompileTime,
MatrixType::MaxColsAtCompileTime > 
ImageResultType
typedef Matrix< int,
MatrixType::RowsAtCompileTime,
1, MatrixType::Options,
MatrixType::MaxRowsAtCompileTime, 1 > 
IntColVectorType
typedef Matrix< int,
1, MatrixType::ColsAtCompileTime,
MatrixType::Options,
1, MatrixType::MaxColsAtCompileTime > 
IntRowVectorType
typedef Matrix< typename
MatrixType::Scalar,
MatrixType::ColsAtCompileTime,
Dynamic, MatrixType::Options,
MatrixType::MaxColsAtCompileTime,
MatrixType::MaxColsAtCompileTime > 
KernelResultType
typedef NumTraits< typename
MatrixType::Scalar >::Real 
RealScalar
typedef Matrix< Scalar,
1, MatrixType::ColsAtCompileTime,
MatrixType::Options,
1, MatrixType::MaxColsAtCompileTime > 
RowVectorType
typedef MatrixType::Scalar Scalar

Public Member Functions

template<typename ImageMatrixType >
void computeImage (ImageMatrixType *result) const
template<typename ResultType >
void computeInverse (ResultType *result) const
template<typename KernelMatrixType >
void computeKernel (KernelMatrixType *result) const
ei_traits< MatrixType >::Scalar determinant () const
int dimensionOfKernel () const
const ImageResultType image () const
MatrixType inverse () const
bool isInjective () const
bool isInvertible () const
bool isSurjective () const
const KernelResultType kernel () const
 LU (const MatrixType &matrix)
const MatrixType & matrixLU () const
const IntColVectorTypepermutationP () const
const IntRowVectorTypepermutationQ () const
int rank () const
template<typename OtherDerived , typename ResultType >
bool solve (const MatrixBase< OtherDerived > &b, ResultType *result) const

Protected Attributes

int m_det_pq
MatrixType m_lu
const MatrixType & m_originalMatrix
IntColVectorType m_p
RealScalar m_precision
IntRowVectorType m_q
int m_rank

Detailed Description

template<typename MatrixType>
class LU< MatrixType >

LU decomposition of a matrix with complete pivoting, and related features.

Parameters:
MatrixType the type of the matrix of which we are computing the LU decomposition

This class represents a LU decomposition of any matrix, with complete pivoting: the matrix A is decomposed as A = PLUQ where L is unit-lower-triangular, U is upper-triangular, and P and Q are permutation matrices. This is a rank-revealing LU decomposition. The eigenvalues (diagonal coefficients) of U are sorted in such a way that any zeros are at the end, so that the rank of A is the index of the first zero on the diagonal of U (with indices starting at 0) if any.

This decomposition provides the generic approach to solving systems of linear equations, computing the rank, invertibility, inverse, kernel, and determinant.

This LU decomposition is very stable and well tested with large matrices. Even exact rank computation works at sizes larger than 1000x1000. However there are use cases where the SVD decomposition is inherently more stable when dealing with numerically damaged input. For example, computing the kernel is more stable with SVD because the SVD can determine which singular values are negligible while LU has to work at the level of matrix coefficients that are less meaningful in this respect.

The data of the LU decomposition can be directly accessed through the methods matrixLU(), permutationP(), permutationQ().

As an exemple, here is how the original matrix can be retrieved:

Output:

See also:
MatrixBase::lu(), MatrixBase::determinant(), MatrixBase::inverse(), MatrixBase::computeInverse()

Definition at line 60 of file LU.h.


Member Typedef Documentation

template<typename MatrixType>
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1, MatrixType::Options, MatrixType::MaxRowsAtCompileTime, 1> LU< MatrixType >::ColVectorType

Definition at line 69 of file LU.h.

template<typename MatrixType>
typedef Matrix<typename MatrixType::Scalar, MatrixType::RowsAtCompileTime, Dynamic, MatrixType::Options, MatrixType::MaxRowsAtCompileTime, MatrixType::MaxColsAtCompileTime > LU< MatrixType >::ImageResultType

Definition at line 93 of file LU.h.

template<typename MatrixType>
typedef Matrix<int, MatrixType::RowsAtCompileTime, 1, MatrixType::Options, MatrixType::MaxRowsAtCompileTime, 1> LU< MatrixType >::IntColVectorType

Definition at line 67 of file LU.h.

template<typename MatrixType>
typedef Matrix<int, 1, MatrixType::ColsAtCompileTime, MatrixType::Options, 1, MatrixType::MaxColsAtCompileTime> LU< MatrixType >::IntRowVectorType

Definition at line 66 of file LU.h.

template<typename MatrixType>
typedef Matrix<typename MatrixType::Scalar, MatrixType::ColsAtCompileTime, Dynamic, MatrixType::Options, MatrixType::MaxColsAtCompileTime, MatrixType::MaxColsAtCompileTime > LU< MatrixType >::KernelResultType

Definition at line 84 of file LU.h.

template<typename MatrixType>
typedef NumTraits<typename MatrixType::Scalar>::Real LU< MatrixType >::RealScalar

Definition at line 65 of file LU.h.

template<typename MatrixType>
typedef Matrix<Scalar, 1, MatrixType::ColsAtCompileTime, MatrixType::Options, 1, MatrixType::MaxColsAtCompileTime> LU< MatrixType >::RowVectorType

Definition at line 68 of file LU.h.

template<typename MatrixType>
typedef MatrixType::Scalar LU< MatrixType >::Scalar

Definition at line 64 of file LU.h.


Member Enumeration Documentation

template<typename MatrixType>
anonymous enum
Enumerator:
MaxSmallDimAtCompileTime 

Definition at line 71 of file LU.h.


Constructor & Destructor Documentation

template<typename MatrixType >
LU< MatrixType >::LU ( const MatrixType &  matrix  )  [inline]

Constructor.

Parameters:
matrix the matrix of which to compute the LU decomposition.

Definition at line 331 of file LU.h.


Member Function Documentation

template<typename MatrixType >
template<typename ImageMatrixType >
void LU< MatrixType >::computeImage ( ImageMatrixType *  result  )  const [inline]

Computes a basis of the image of the matrix, also called the column-space or range of he matrix.

Note:
Calling this method on the zero matrix will make an assertion fail.
Parameters:
result a pointer to the matrix in which to store the image. The columns of this matrix will be set to form a basis of the image (it will be resized if necessary).

Example:

Output:

See also:
image(), computeKernel(), kernel()

Definition at line 457 of file LU.h.

template<typename MatrixType>
template<typename ResultType >
void LU< MatrixType >::computeInverse ( ResultType *  result  )  const [inline]

Computes the inverse of the matrix of which *this is the LU decomposition.

Parameters:
result a pointer to the matrix into which to store the inverse. Resized if needed.
Note:
If this matrix is not invertible, *result is left with undefined coefficients. Use isInvertible() to first determine whether this matrix is invertible.
See also:
MatrixBase::computeInverse(), inverse()

Definition at line 301 of file LU.h.

template<typename MatrixType >
template<typename KernelMatrixType >
void LU< MatrixType >::computeKernel ( KernelMatrixType *  result  )  const [inline]

Computes a basis of the kernel of the matrix, also called the null-space of the matrix.

Note:
This method is only allowed on non-invertible matrices, as determined by isInvertible(). Calling it on an invertible matrix will make an assertion fail.
Parameters:
result a pointer to the matrix in which to store the kernel. The columns of this matrix will be set to form a basis of the kernel (it will be resized if necessary).

Example:

Output:

See also:
kernel(), computeImage(), image()

Definition at line 411 of file LU.h.

template<typename MatrixType >
ei_traits< MatrixType >::Scalar LU< MatrixType >::determinant (  )  const [inline]
Returns:
the determinant of the matrix of which *this is the LU decomposition. It has only linear complexity (that is, O(n) where n is the dimension of the square matrix) as the LU decomposition has already been computed.
Note:
This is only for square matrices.
For fixed-size matrices of size up to 4, MatrixBase::determinant() offers optimized paths.
Warning:
a determinant can be very big or small, so for matrices of large enough dimension, there is a risk of overflow/underflow.
See also:
MatrixBase::determinant()

Definition at line 404 of file LU.h.

template<typename MatrixType>
int LU< MatrixType >::dimensionOfKernel (  )  const [inline]
Returns:
the dimension of the kernel of the matrix of which *this is the LU decomposition.
Note:
Since the rank is computed at the time of the construction of the LU decomposition, this method almost does not perform any further computation.

Definition at line 254 of file LU.h.

template<typename MatrixType >
const LU< MatrixType >::ImageResultType LU< MatrixType >::image (  )  const [inline]
Returns:
the image of the matrix, also called its column-space. The columns of the returned matrix will form a basis of the kernel.
Note:
: Calling this method on the zero matrix will make an assertion fail.
: this method returns a matrix by value, which induces some inefficiency. If you prefer to avoid this overhead, use computeImage() instead.

Example:

Output:

See also:
computeImage(), kernel()

Definition at line 467 of file LU.h.

template<typename MatrixType>
MatrixType LU< MatrixType >::inverse ( void   )  const [inline]
Returns:
the inverse of the matrix of which *this is the LU decomposition.
Note:
If this matrix is not invertible, the returned matrix has undefined coefficients. Use isInvertible() to first determine whether this matrix is invertible.
See also:
computeInverse(), MatrixBase::inverse()

Definition at line 313 of file LU.h.

template<typename MatrixType>
bool LU< MatrixType >::isInjective (  )  const [inline]
Returns:
true if the matrix of which *this is the LU decomposition represents an injective linear map, i.e. has trivial kernel; false otherwise.
Note:
Since the rank is computed at the time of the construction of the LU decomposition, this method almost does not perform any further computation.

Definition at line 265 of file LU.h.

template<typename MatrixType>
bool LU< MatrixType >::isInvertible (  )  const [inline]
Returns:
true if the matrix of which *this is the LU decomposition is invertible.
Note:
Since the rank is computed at the time of the construction of the LU decomposition, this method almost does not perform any further computation.

Definition at line 286 of file LU.h.

template<typename MatrixType>
bool LU< MatrixType >::isSurjective (  )  const [inline]
Returns:
true if the matrix of which *this is the LU decomposition represents a surjective linear map; false otherwise.
Note:
Since the rank is computed at the time of the construction of the LU decomposition, this method almost does not perform any further computation.

Definition at line 276 of file LU.h.

template<typename MatrixType >
const LU< MatrixType >::KernelResultType LU< MatrixType >::kernel (  )  const [inline]
Returns:
the kernel of the matrix, also called its null-space. The columns of the returned matrix will form a basis of the kernel.
Note:
: this method is only allowed on non-invertible matrices, as determined by isInvertible(). Calling it on an invertible matrix will make an assertion fail.
: this method returns a matrix by value, which induces some inefficiency. If you prefer to avoid this overhead, use computeKernel() instead.

Example:

Output:

See also:
computeKernel(), image()

Definition at line 448 of file LU.h.

template<typename MatrixType>
const MatrixType& LU< MatrixType >::matrixLU (  )  const [inline]
Returns:
the LU decomposition matrix: the upper-triangular part is U, the unit-lower-triangular part is L (at least for square matrices; in the non-square case, special care is needed, see the documentation of class LU).
See also:
matrixL(), matrixU()

Definition at line 107 of file LU.h.

template<typename MatrixType>
const IntColVectorType& LU< MatrixType >::permutationP (  )  const [inline]
Returns:
a vector of integers, whose size is the number of rows of the matrix being decomposed, representing the P permutation i.e. the permutation of the rows. For its precise meaning, see the examples given in the documentation of class LU.
See also:
permutationQ()

Definition at line 118 of file LU.h.

template<typename MatrixType>
const IntRowVectorType& LU< MatrixType >::permutationQ (  )  const [inline]
Returns:
a vector of integers, whose size is the number of columns of the matrix being decomposed, representing the Q permutation i.e. the permutation of the columns. For its precise meaning, see the examples given in the documentation of class LU.
See also:
permutationP()

Definition at line 129 of file LU.h.

template<typename MatrixType>
int LU< MatrixType >::rank (  )  const [inline]
Returns:
the rank of the matrix of which *this is the LU decomposition.
Note:
This is computed at the time of the construction of the LU decomposition. This method does not perform any further computation.

Definition at line 244 of file LU.h.

template<typename MatrixType >
template<typename OtherDerived , typename ResultType >
bool LU< MatrixType >::solve ( const MatrixBase< OtherDerived > &  b,
ResultType *  result 
) const [inline]

This method finds a solution x to the equation Ax=b, where A is the matrix of which *this is the LU decomposition, if any exists.

Parameters:
b the right-hand-side of the equation to solve. Can be a vector or a matrix, the only requirement in order for the equation to make sense is that b.rows()==A.rows(), where A is the matrix of which *this is the LU decomposition.
result a pointer to the vector or matrix in which to store the solution, if any exists. Resized if necessary, so that result->rows()==A.cols() and result->cols()==b.cols(). If no solution exists, *result is left with undefined coefficients.
Returns:
true if any solution exists, false if no solution exists.
Note:
If there exist more than one solution, this method will arbitrarily choose one. If you need a complete analysis of the space of solutions, take the one solution obtained by this method and add to it elements of the kernel, as determined by kernel().

Example:

Output:

See also:
MatrixBase::solveTriangular(), kernel(), computeKernel(), inverse(), computeInverse()

Definition at line 476 of file LU.h.


Member Data Documentation

template<typename MatrixType>
int LU< MatrixType >::m_det_pq [protected]

Definition at line 325 of file LU.h.

template<typename MatrixType>
MatrixType LU< MatrixType >::m_lu [protected]

Definition at line 322 of file LU.h.

template<typename MatrixType>
const MatrixType& LU< MatrixType >::m_originalMatrix [protected]

Definition at line 321 of file LU.h.

template<typename MatrixType>
IntColVectorType LU< MatrixType >::m_p [protected]

Definition at line 323 of file LU.h.

template<typename MatrixType>
RealScalar LU< MatrixType >::m_precision [protected]

Definition at line 327 of file LU.h.

template<typename MatrixType>
IntRowVectorType LU< MatrixType >::m_q [protected]

Definition at line 324 of file LU.h.

template<typename MatrixType>
int LU< MatrixType >::m_rank [protected]

Definition at line 326 of file LU.h.


The documentation for this class was generated from the following file:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines


vcglib
Author(s): Christian Bersch
autogenerated on Fri Jan 11 09:22:05 2013