Public Types | Public Member Functions | Protected Attributes
FullPivLU< _MatrixType > Class Template Reference

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

#include <FullPivLU.h>

List of all members.

Public Types

enum  {
  RowsAtCompileTime = MatrixType::RowsAtCompileTime, ColsAtCompileTime = MatrixType::ColsAtCompileTime, Options = MatrixType::Options, MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
}
typedef MatrixType::Index Index
typedef
internal::plain_col_type
< MatrixType, Index >::type 
IntColVectorType
typedef
internal::plain_row_type
< MatrixType, Index >::type 
IntRowVectorType
typedef _MatrixType MatrixType
typedef PermutationMatrix
< RowsAtCompileTime,
MaxRowsAtCompileTime
PermutationPType
typedef PermutationMatrix
< ColsAtCompileTime,
MaxColsAtCompileTime
PermutationQType
typedef NumTraits< typename
MatrixType::Scalar >::Real 
RealScalar
typedef MatrixType::Scalar Scalar
typedef internal::traits
< MatrixType >::StorageKind 
StorageKind

Public Member Functions

Index cols () const
FullPivLUcompute (const MatrixType &matrix)
internal::traits< MatrixType >
::Scalar 
determinant () const
Index dimensionOfKernel () const
 FullPivLU ()
 Default Constructor.
 FullPivLU (Index rows, Index cols)
 Default Constructor with memory preallocation.
 FullPivLU (const MatrixType &matrix)
const internal::image_retval
< FullPivLU
image (const MatrixType &originalMatrix) const
const internal::solve_retval
< FullPivLU, typename
MatrixType::IdentityReturnType > 
inverse () const
bool isInjective () const
bool isInvertible () const
bool isSurjective () const
const internal::kernel_retval
< FullPivLU
kernel () const
const MatrixTypematrixLU () const
RealScalar maxPivot () const
Index nonzeroPivots () const
const PermutationPTypepermutationP () const
const PermutationQTypepermutationQ () const
Index rank () const
MatrixType reconstructedMatrix () const
Index rows () const
FullPivLUsetThreshold (const RealScalar &threshold)
FullPivLUsetThreshold (Default_t)
template<typename Rhs >
const internal::solve_retval
< FullPivLU, Rhs > 
solve (const MatrixBase< Rhs > &b) const
RealScalar threshold () const

Protected Attributes

IntRowVectorType m_colsTranspositions
Index m_det_pq
bool m_isInitialized
MatrixType m_lu
RealScalar m_maxpivot
Index m_nonzero_pivots
PermutationPType m_p
RealScalar m_prescribedThreshold
PermutationQType m_q
IntColVectorType m_rowsTranspositions
bool m_usePrescribedThreshold

Detailed Description

template<typename _MatrixType>
class FullPivLU< _MatrixType >

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

Parameters:
MatrixTypethe 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.

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. However there are use cases where the SVD decomposition is inherently more stable and/or flexible. For example, when computing the kernel of a matrix, working with the SVD allows to select the smallest singular values of the matrix, something that the LU decomposition doesn't see.

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:

typedef Matrix<double, 5, 3> Matrix5x3;
typedef Matrix<double, 5, 5> Matrix5x5;
Matrix5x3 m = Matrix5x3::Random();
cout << "Here is the matrix m:" << endl << m << endl;
Eigen::FullPivLU<Matrix5x3> lu(m);
cout << "Here is, up to permutations, its LU decomposition matrix:"
     << endl << lu.matrixLU() << endl;
cout << "Here is the L part:" << endl;
Matrix5x5 l = Matrix5x5::Identity();
l.block<5,3>(0,0).triangularView<StrictlyLower>() = lu.matrixLU();
cout << l << endl;
cout << "Here is the U part:" << endl;
Matrix5x3 u = lu.matrixLU().triangularView<Upper>();
cout << u << endl;
cout << "Let us now reconstruct the original matrix m:" << endl;
cout << lu.permutationP().inverse() * l * u * lu.permutationQ().inverse() << endl;

Output:

See also:
MatrixBase::fullPivLu(), MatrixBase::determinant(), MatrixBase::inverse()

Definition at line 58 of file FullPivLU.h.


Member Typedef Documentation

template<typename _MatrixType>
typedef MatrixType::Index FullPivLU< _MatrixType >::Index

Definition at line 72 of file FullPivLU.h.

template<typename _MatrixType>
typedef internal::plain_col_type<MatrixType, Index>::type FullPivLU< _MatrixType >::IntColVectorType

Reimplemented in LU< MatrixType >.

Definition at line 74 of file FullPivLU.h.

template<typename _MatrixType>
typedef internal::plain_row_type<MatrixType, Index>::type FullPivLU< _MatrixType >::IntRowVectorType

Reimplemented in LU< MatrixType >.

Definition at line 73 of file FullPivLU.h.

template<typename _MatrixType>
typedef _MatrixType FullPivLU< _MatrixType >::MatrixType

Definition at line 61 of file FullPivLU.h.

template<typename _MatrixType>
typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime> FullPivLU< _MatrixType >::PermutationPType

Definition at line 76 of file FullPivLU.h.

template<typename _MatrixType>
typedef PermutationMatrix<ColsAtCompileTime, MaxColsAtCompileTime> FullPivLU< _MatrixType >::PermutationQType

Definition at line 75 of file FullPivLU.h.

template<typename _MatrixType>
typedef NumTraits<typename MatrixType::Scalar>::Real FullPivLU< _MatrixType >::RealScalar

Reimplemented in LU< MatrixType >.

Definition at line 70 of file FullPivLU.h.

template<typename _MatrixType>
typedef MatrixType::Scalar FullPivLU< _MatrixType >::Scalar

Reimplemented in LU< MatrixType >.

Definition at line 69 of file FullPivLU.h.

template<typename _MatrixType>
typedef internal::traits<MatrixType>::StorageKind FullPivLU< _MatrixType >::StorageKind

Definition at line 71 of file FullPivLU.h.


Member Enumeration Documentation

template<typename _MatrixType>
anonymous enum
Enumerator:
RowsAtCompileTime 
ColsAtCompileTime 
Options 
MaxRowsAtCompileTime 
MaxColsAtCompileTime 

Definition at line 62 of file FullPivLU.h.


Constructor & Destructor Documentation

template<typename _MatrixType>
FullPivLU< _MatrixType >::FullPivLU ( )

Default Constructor.

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

template<typename _MatrixType>
FullPivLU< _MatrixType >::FullPivLU ( Index  rows,
Index  cols 
)

Default Constructor with memory preallocation.

Like the default constructor but with preallocation of the internal data according to the specified problem size.

See also:
FullPivLU()
template<typename _MatrixType>
FullPivLU< _MatrixType >::FullPivLU ( const MatrixType matrix)

Constructor.

Parameters:
matrixthe matrix of which to compute the LU decomposition. It is required to be nonzero.

Member Function Documentation

template<typename _MatrixType>
Index FullPivLU< _MatrixType >::cols ( void  ) const [inline]

Definition at line 384 of file FullPivLU.h.

template<typename _MatrixType>
FullPivLU& FullPivLU< _MatrixType >::compute ( const MatrixType matrix)

Computes the LU decomposition of the given matrix.

Parameters:
matrixthe matrix of which to compute the LU decomposition. It is required to be nonzero.
Returns:
a reference to *this
template<typename _MatrixType>
internal::traits<MatrixType>::Scalar FullPivLU< _MatrixType >::determinant ( ) const
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()
template<typename _MatrixType>
Index FullPivLU< _MatrixType >::dimensionOfKernel ( ) const [inline]
Returns:
the dimension of the kernel of the matrix of which *this is the LU decomposition.
Note:
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).

Definition at line 322 of file FullPivLU.h.

template<typename _MatrixType>
const internal::image_retval<FullPivLU> FullPivLU< _MatrixType >::image ( const MatrixType originalMatrix) 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.
Parameters:
originalMatrixthe original matrix, of which *this is the LU decomposition. The reason why it is needed to pass it here, is that this allows a large optimization, as otherwise this method would need to reconstruct it from the LU decomposition.
Note:
If the image has dimension zero, then the returned matrix is a column-vector filled with zeros.
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).

Example:

Matrix3d m;
m << 1,1,0,
     1,3,2,
     0,1,1;
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Notice that the middle column is the sum of the two others, so the "
     << "columns are linearly dependent." << endl;
cout << "Here is a matrix whose columns have the same span but are linearly independent:"
     << endl << m.fullPivLu().image(m) << endl;

Output:

See also:
kernel()

Definition at line 200 of file FullPivLU.h.

template<typename _MatrixType>
const internal::solve_retval<FullPivLU,typename MatrixType::IdentityReturnType> FullPivLU< _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:
MatrixBase::inverse()

Definition at line 373 of file FullPivLU.h.

template<typename _MatrixType>
bool FullPivLU< _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:
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).

Definition at line 335 of file FullPivLU.h.

template<typename _MatrixType>
bool FullPivLU< _MatrixType >::isInvertible ( ) const [inline]
Returns:
true if the matrix of which *this is the LU decomposition is invertible.
Note:
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).

Definition at line 360 of file FullPivLU.h.

template<typename _MatrixType>
bool FullPivLU< _MatrixType >::isSurjective ( ) const [inline]
Returns:
true if the matrix of which *this is the LU decomposition represents a surjective linear map; false otherwise.
Note:
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).

Definition at line 348 of file FullPivLU.h.

template<typename _MatrixType>
const internal::kernel_retval<FullPivLU> FullPivLU< _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:
If the kernel has dimension zero, then the returned matrix is a column-vector filled with zeros.
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).

Example:

MatrixXf m = MatrixXf::Random(3,5);
cout << "Here is the matrix m:" << endl << m << endl;
MatrixXf ker = m.fullPivLu().kernel();
cout << "Here is a matrix whose columns form a basis of the kernel of m:"
     << endl << ker << endl;
cout << "By definition of the kernel, m*ker is zero:"
     << endl << m*ker << endl;

Output:

See also:
image()

Definition at line 174 of file FullPivLU.h.

template<typename _MatrixType>
const MatrixType& FullPivLU< _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 FullPivLU).
See also:
matrixL(), matrixU()

Definition at line 116 of file FullPivLU.h.

template<typename _MatrixType>
RealScalar FullPivLU< _MatrixType >::maxPivot ( ) const [inline]
Returns:
the absolute value of the biggest pivot, i.e. the biggest diagonal coefficient of U.

Definition at line 138 of file FullPivLU.h.

template<typename _MatrixType>
Index FullPivLU< _MatrixType >::nonzeroPivots ( ) const [inline]
Returns:
the number of nonzero pivots in the LU decomposition. Here nonzero is meant in the exact sense, not in a fuzzy sense. So that notion isn't really intrinsically interesting, but it is still useful when implementing algorithms.
See also:
rank()

Definition at line 129 of file FullPivLU.h.

template<typename _MatrixType>
const PermutationPType& FullPivLU< _MatrixType >::permutationP ( ) const [inline]
Returns:
the permutation matrix P
See also:
permutationQ()

Definition at line 144 of file FullPivLU.h.

template<typename _MatrixType>
const PermutationQType& FullPivLU< _MatrixType >::permutationQ ( ) const [inline]
Returns:
the permutation matrix Q
See also:
permutationP()

Definition at line 154 of file FullPivLU.h.

template<typename _MatrixType>
Index FullPivLU< _MatrixType >::rank ( ) const [inline]
Returns:
the rank of the matrix of which *this is the LU decomposition.
Note:
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).

Definition at line 306 of file FullPivLU.h.

template<typename _MatrixType>
MatrixType FullPivLU< _MatrixType >::reconstructedMatrix ( ) const
template<typename _MatrixType>
Index FullPivLU< _MatrixType >::rows ( void  ) const [inline]

Definition at line 383 of file FullPivLU.h.

template<typename _MatrixType>
FullPivLU& FullPivLU< _MatrixType >::setThreshold ( const RealScalar threshold) [inline]

Allows to prescribe a threshold to be used by certain methods, such as rank(), who need to determine when pivots are to be considered nonzero. This is not used for the LU decomposition itself.

When it needs to get the threshold value, Eigen calls threshold(). By default, this uses a formula to automatically determine a reasonable threshold. Once you have called the present method setThreshold(const RealScalar&), your value is used instead.

Parameters:
thresholdThe new value to use as the threshold.

A pivot will be considered nonzero if its absolute value is strictly greater than $ \vert pivot \vert \leqslant threshold \times \vert maxpivot \vert $ where maxpivot is the biggest pivot.

If you want to come back to the default behavior, call setThreshold(Default_t)

Definition at line 267 of file FullPivLU.h.

template<typename _MatrixType>
FullPivLU& FullPivLU< _MatrixType >::setThreshold ( Default_t  ) [inline]

Allows to come back to the default behavior, letting Eigen use its default formula for determining the threshold.

You should pass the special object Eigen::Default as parameter here.

 lu.setThreshold(Eigen::Default); 

See the documentation of setThreshold(const RealScalar&).

Definition at line 282 of file FullPivLU.h.

template<typename _MatrixType>
template<typename Rhs >
const internal::solve_retval<FullPivLU, Rhs> FullPivLU< _MatrixType >::solve ( const MatrixBase< Rhs > &  b) const [inline]
Returns:
a solution x to the equation Ax=b, where A is the matrix of which *this is the LU decomposition.
Parameters:
bthe 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.
Returns:
a solution.

Example:

Matrix<float,2,3> m = Matrix<float,2,3>::Random();
Matrix2f y = Matrix2f::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is the matrix y:" << endl << y << endl;
Matrix<float,3,2> x = m.fullPivLu().solve(y);
if((m*x).isApprox(y))
{
  cout << "Here is a solution x to the equation mx=y:" << endl << x << endl;
}
else
  cout << "The equation mx=y does not have any solution." << endl;

Output:

See also:
TriangularView::solve(), kernel(), inverse()

Definition at line 227 of file FullPivLU.h.

template<typename _MatrixType>
RealScalar FullPivLU< _MatrixType >::threshold ( ) const [inline]

Returns the threshold that will be used by certain methods such as rank().

See the documentation of setThreshold(const RealScalar&).

Definition at line 291 of file FullPivLU.h.


Member Data Documentation

template<typename _MatrixType>
IntRowVectorType FullPivLU< _MatrixType >::m_colsTranspositions [protected]

Definition at line 391 of file FullPivLU.h.

template<typename _MatrixType>
Index FullPivLU< _MatrixType >::m_det_pq [protected]

Definition at line 392 of file FullPivLU.h.

template<typename _MatrixType>
bool FullPivLU< _MatrixType >::m_isInitialized [protected]

Definition at line 394 of file FullPivLU.h.

template<typename _MatrixType>
MatrixType FullPivLU< _MatrixType >::m_lu [protected]

Definition at line 387 of file FullPivLU.h.

template<typename _MatrixType>
RealScalar FullPivLU< _MatrixType >::m_maxpivot [protected]

Definition at line 393 of file FullPivLU.h.

template<typename _MatrixType>
Index FullPivLU< _MatrixType >::m_nonzero_pivots [protected]

Definition at line 392 of file FullPivLU.h.

template<typename _MatrixType>
PermutationPType FullPivLU< _MatrixType >::m_p [protected]

Definition at line 388 of file FullPivLU.h.

template<typename _MatrixType>
RealScalar FullPivLU< _MatrixType >::m_prescribedThreshold [protected]

Definition at line 393 of file FullPivLU.h.

template<typename _MatrixType>
PermutationQType FullPivLU< _MatrixType >::m_q [protected]

Definition at line 389 of file FullPivLU.h.

template<typename _MatrixType>
IntColVectorType FullPivLU< _MatrixType >::m_rowsTranspositions [protected]

Definition at line 390 of file FullPivLU.h.

template<typename _MatrixType>
bool FullPivLU< _MatrixType >::m_usePrescribedThreshold [protected]

Definition at line 394 of file FullPivLU.h.


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


libicr
Author(s): Robert Krug
autogenerated on Mon Jan 6 2014 11:34:15