Eigen  3.4.90 (git rev 67eeba6e720c5745abc77ae6c92ce0a44aa7b7ae)
Eigen::FullPivLU< MatrixType_ > Class Template Reference

Detailed Description

template<typename MatrixType_>
class Eigen::FullPivLU< MatrixType_ >

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

Template 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 = P^{-1} L U Q^{-1} \) 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 example, 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;
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;
LU decomposition of a matrix with complete pivoting, and related features.
Definition: FullPivLU.h:64
@ Upper
Definition: Constants.h:213

Output:

Here is the matrix m:
   0.68  -0.605 -0.0452
 -0.211   -0.33   0.258
  0.566   0.536   -0.27
  0.597  -0.444  0.0268
  0.823   0.108   0.904
Here is, up to permutations, its LU decomposition matrix:
 0.904  0.823  0.108
-0.299  0.812  0.569
 -0.05  0.888   -1.1
0.0296  0.705  0.768
 0.285 -0.549 0.0436
Here is the L part:
     1      0      0      0      0
-0.299      1      0      0      0
 -0.05  0.888      1      0      0
0.0296  0.705  0.768      1      0
 0.285 -0.549 0.0436      0      1
Here is the U part:
0.904 0.823 0.108
    0 0.812 0.569
    0     0  -1.1
    0     0     0
    0     0     0
Let us now reconstruct the original matrix m:
   0.68  -0.605 -0.0452
 -0.211   -0.33   0.258
  0.566   0.536   -0.27
  0.597  -0.444  0.0268
  0.823   0.108   0.904

This class supports the inplace decomposition mechanism.

See also
MatrixBase::fullPivLu(), MatrixBase::determinant(), MatrixBase::inverse()
+ Inheritance diagram for Eigen::FullPivLU< MatrixType_ >:

Public Member Functions

template<typename InputType >
FullPivLUcompute (const EigenBase< InputType > &matrix)
 
internal::traits< MatrixType >::Scalar determinant () const
 
Index dimensionOfKernel () const
 
 FullPivLU ()
 Default Constructor. More...
 
template<typename InputType >
 FullPivLU (const EigenBase< InputType > &matrix)
 
template<typename InputType >
 FullPivLU (EigenBase< InputType > &matrix)
 Constructs a LU factorization from a given matrix. More...
 
 FullPivLU (Index rows, Index cols)
 Default Constructor with memory preallocation. More...
 
const internal::image_retval< FullPivLUimage (const MatrixType &originalMatrix) const
 
const Inverse< FullPivLUinverse () const
 
bool isInjective () const
 
bool isInvertible () const
 
bool isSurjective () const
 
const internal::kernel_retval< FullPivLUkernel () const
 
const MatrixType & matrixLU () const
 
RealScalar maxPivot () const
 
Index nonzeroPivots () const
 
const PermutationPTypepermutationP () const
 
const PermutationQTypepermutationQ () const
 
Index rank () const
 
RealScalar rcond () const
 
MatrixType reconstructedMatrix () const
 
FullPivLUsetThreshold (const RealScalar &threshold)
 
FullPivLUsetThreshold (Default_t)
 
template<typename Rhs >
const Solve< FullPivLU, Rhs > solve (const MatrixBase< Rhs > &b) const
 
RealScalar threshold () const
 
- Public Member Functions inherited from Eigen::SolverBase< FullPivLU< MatrixType_ > >
const AdjointReturnType adjoint () const
 
FullPivLU< MatrixType_ > & derived ()
 
const FullPivLU< MatrixType_ > & derived () const
 
const Solve< FullPivLU< MatrixType_ >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
 SolverBase ()
 
const ConstTransposeReturnType transpose () const
 
- Public Member Functions inherited from Eigen::EigenBase< Derived >
EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
Derived & derived ()
 
const Derived & derived () const
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
EIGEN_CONSTEXPR Index size () const EIGEN_NOEXCEPT
 

Additional Inherited Members

- Public Types inherited from Eigen::EigenBase< Derived >
typedef Eigen::Index Index
 The interface type of indices. More...
 

Constructor & Destructor Documentation

◆ FullPivLU() [1/4]

template<typename MatrixType >
Eigen::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&).

◆ FullPivLU() [2/4]

template<typename MatrixType >
Eigen::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()

◆ FullPivLU() [3/4]

template<typename MatrixType >
template<typename InputType >
Eigen::FullPivLU< MatrixType >::FullPivLU ( const EigenBase< InputType > &  matrix)
explicit

Constructor.

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

◆ FullPivLU() [4/4]

template<typename MatrixType >
template<typename InputType >
Eigen::FullPivLU< MatrixType >::FullPivLU ( EigenBase< InputType > &  matrix)
explicit

Constructs a LU factorization from a given matrix.

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

See also
FullPivLU(const EigenBase&)

Member Function Documentation

◆ compute()

template<typename MatrixType_ >
template<typename InputType >
FullPivLU& Eigen::FullPivLU< MatrixType_ >::compute ( const EigenBase< InputType > &  matrix)
inline

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

◆ determinant()

template<typename MatrixType >
internal::traits< MatrixType >::Scalar Eigen::FullPivLU< MatrixType >::determinant
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()

◆ dimensionOfKernel()

template<typename MatrixType_ >
Index Eigen::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&).

◆ image()

template<typename MatrixType_ >
const internal::image_retval<FullPivLU> Eigen::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 image (column-space).
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:

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;
Matrix< double, 3, 3 > Matrix3d
3×3 matrix of type double.
Definition: Matrix.h:501

Output:

Here is the matrix m:
1 1 0
1 3 2
0 1 1
Notice that the middle column is the sum of the two others, so the columns are linearly dependent.
Here is a matrix whose columns have the same span but are linearly independent:
1 1
3 1
1 0
See also
kernel()

◆ inverse()

template<typename MatrixType_ >
const Inverse<FullPivLU> Eigen::FullPivLU< MatrixType_ >::inverse ( ) 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()

◆ isInjective()

template<typename MatrixType_ >
bool Eigen::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&).

◆ isInvertible()

template<typename MatrixType_ >
bool Eigen::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&).

◆ isSurjective()

template<typename MatrixType_ >
bool Eigen::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&).

◆ kernel()

template<typename MatrixType_ >
const internal::kernel_retval<FullPivLU> Eigen::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:

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;
static const RandomReturnType Random()
Definition: Random.h:114
Matrix< float, Dynamic, Dynamic > MatrixXf
Dynamic×Dynamic matrix of type float.
Definition: Matrix.h:500

Output:

Here is the matrix m:
   0.68   0.597   -0.33   0.108   -0.27
 -0.211   0.823   0.536 -0.0452  0.0268
  0.566  -0.605  -0.444   0.258   0.904
Here is a matrix whose columns form a basis of the kernel of m:
 -0.219   0.763
0.00335  -0.447
      0       1
      1       0
 -0.145  -0.285
By definition of the kernel, m*ker is zero:
 7.45e-09  1.49e-08
-1.86e-09 -4.05e-08
        0 -2.98e-08
See also
image()

◆ matrixLU()

template<typename MatrixType_ >
const MatrixType& Eigen::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()

◆ maxPivot()

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

◆ nonzeroPivots()

template<typename MatrixType_ >
Index Eigen::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()

◆ permutationP()

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

◆ permutationQ()

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

◆ rank()

template<typename MatrixType_ >
Index Eigen::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&).

◆ rcond()

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

◆ reconstructedMatrix()

template<typename MatrixType >
MatrixType Eigen::FullPivLU< MatrixType >::reconstructedMatrix
Returns
the matrix represented by the decomposition, i.e., it returns the product: \( P^{-1} L U Q^{-1} \). This function is provided for debug purposes.

◆ setThreshold() [1/2]

template<typename MatrixType_ >
FullPivLU& Eigen::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)

◆ setThreshold() [2/2]

template<typename MatrixType_ >
FullPivLU& Eigen::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&).

◆ solve()

template<typename MatrixType_ >
template<typename Rhs >
const Solve<FullPivLU, Rhs> Eigen::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.

This method just tries to find as good a solution as possible. If you want to check whether a solution exists or if it is accurate, just call this function to get a result and then compute the error of this result, or use MatrixBase::isApprox() directly, for instance like this:

bool a_solution_exists = (A*result).isApprox(b, precision);

This method avoids dividing by zero, so that the non-existence of a solution doesn't by itself mean that you'll get inf or nan values.

If there exists 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:

Matrix<float,2,3> m = Matrix<float,2,3>::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;
Matrix< float, 2, 2 > Matrix2f
2×2 matrix of type float.
Definition: Matrix.h:500

Output:

Here is the matrix m:
  0.68  0.566  0.823
-0.211  0.597 -0.605
Here is the matrix y:
 -0.33 -0.444
 0.536  0.108
Here is a solution x to the equation mx=y:
     0      0
 0.291 -0.216
  -0.6 -0.391
See also
TriangularView::solve(), kernel(), inverse()

◆ threshold()

template<typename MatrixType_ >
RealScalar Eigen::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&).


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