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

Detailed Description

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

LU decomposition of a matrix with partial 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 a square invertible matrix, with partial pivoting: the matrix A is decomposed as A = PLU where L is unit-lower-triangular, U is upper-triangular, and P is a permutation matrix.

Typically, partial pivoting LU decomposition is only considered numerically stable for square invertible matrices. Thus LAPACK's dgesv and dgesvx require the matrix to be square and invertible. The present class does the same. It will assert that the matrix is square, but it won't (actually it can't) check that the matrix is invertible: it is your task to check that you only use this decomposition on invertible matrices.

The guaranteed safe alternative, working for all matrices, is the full pivoting LU decomposition, provided by class FullPivLU.

This is not a rank-revealing LU decomposition. Many features are intentionally absent from this class, such as rank computation. If you need these features, use class FullPivLU.

This LU decomposition is suitable to invert invertible matrices. It is what MatrixBase::inverse() uses in the general case. On the other hand, it is not suitable to determine whether a given matrix is invertible.

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

This class supports the inplace decomposition mechanism.

See also
MatrixBase::partialPivLu(), MatrixBase::determinant(), MatrixBase::inverse(), MatrixBase::computeInverse(), class FullPivLU
+ Inheritance diagram for Eigen::PartialPivLU< MatrixType_ >:

Public Member Functions

Scalar determinant () const
 
const Inverse< PartialPivLUinverse () const
 
const MatrixType & matrixLU () const
 
 PartialPivLU ()
 Default Constructor. More...
 
template<typename InputType >
 PartialPivLU (const EigenBase< InputType > &matrix)
 
template<typename InputType >
 PartialPivLU (EigenBase< InputType > &matrix)
 
 PartialPivLU (Index size)
 Default Constructor with memory preallocation. More...
 
const PermutationTypepermutationP () const
 
RealScalar rcond () const
 
MatrixType reconstructedMatrix () const
 
template<typename Rhs >
const Solve< PartialPivLU, Rhs > solve (const MatrixBase< Rhs > &b) const
 
- Public Member Functions inherited from Eigen::SolverBase< PartialPivLU< MatrixType_ > >
const AdjointReturnType adjoint () const
 
PartialPivLU< MatrixType_ > & derived ()
 
const PartialPivLU< MatrixType_ > & derived () const
 
const Solve< PartialPivLU< 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

◆ PartialPivLU() [1/4]

template<typename MatrixType >
Eigen::PartialPivLU< MatrixType >::PartialPivLU

Default Constructor.

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

◆ PartialPivLU() [2/4]

template<typename MatrixType >
Eigen::PartialPivLU< MatrixType >::PartialPivLU ( Index  size)
explicit

Default Constructor with memory preallocation.

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

See also
PartialPivLU()

◆ PartialPivLU() [3/4]

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

Constructor.

Parameters
matrixthe matrix of which to compute the LU decomposition.
Warning
The matrix should have full rank (e.g. if it's square, it should be invertible). If you need to deal with non-full rank, use class FullPivLU instead.

◆ PartialPivLU() [4/4]

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

Constructor for inplace decomposition

Parameters
matrixthe matrix of which to compute the LU decomposition.
Warning
The matrix should have full rank (e.g. if it's square, it should be invertible). If you need to deal with non-full rank, use class FullPivLU instead.

Member Function Documentation

◆ determinant()

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

◆ inverse()

template<typename MatrixType_ >
const Inverse<PartialPivLU> Eigen::PartialPivLU< MatrixType_ >::inverse ( ) const
inline
Returns
the inverse of the matrix of which *this is the LU decomposition.
Warning
The matrix being decomposed here is assumed to be invertible. If you need to check for invertibility, use class FullPivLU instead.
See also
MatrixBase::inverse(), LU::inverse()

◆ matrixLU()

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

◆ permutationP()

template<typename MatrixType_ >
const PermutationType& Eigen::PartialPivLU< MatrixType_ >::permutationP ( ) const
inline
Returns
the permutation matrix P.

◆ rcond()

template<typename MatrixType_ >
RealScalar Eigen::PartialPivLU< 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::PartialPivLU< MatrixType >::reconstructedMatrix
Returns
the matrix represented by the decomposition, i.e., it returns the product: P^{-1} L U. This function is provided for debug purpose.

◆ solve()

template<typename MatrixType_ >
template<typename Rhs >
const Solve<PartialPivLU, Rhs> Eigen::PartialPivLU< MatrixType_ >::solve ( const MatrixBase< Rhs > &  b) const
inline

This method returns the 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
the solution.

Example:

cout << "Here is the invertible matrix A:" << endl << A << endl;
cout << "Here is the matrix B:" << endl << B << endl;
MatrixXd X = A.lu().solve(B);
cout << "Here is the (unique) solution X to the equation AX=B:" << endl << X << endl;
cout << "Relative error: " << (A*X-B).norm() / B.norm() << endl;
static const RandomReturnType Random()
Definition: Random.h:114
Matrix< double, Dynamic, Dynamic > MatrixXd
Dynamic×Dynamic matrix of type double.
Definition: Matrix.h:501

Output:

Here is the invertible matrix A:
  0.68  0.597  -0.33
-0.211  0.823  0.536
 0.566 -0.605 -0.444
Here is the matrix B:
  0.108   -0.27
-0.0452  0.0268
  0.258   0.904
Here is the (unique) solution X to the equation AX=B:
 0.609   2.68
-0.231  -1.57
  0.51   3.51
Relative error: 3.28e-16

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

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

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