template<typename _MatrixType>
PartialPivLU class
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 |
Contents
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::
The data of the LU decomposition can be directly accessed through the methods matrixLU(), permutationP().
This class supports the inplace decomposition mechanism.
Base classes
-
template<typename Derived>class SolverBase
- A base class for matrix decomposition and solvers.
Public types
- enum (anonymous) { MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime }
- using Base = SolverBase<PartialPivLU>
- using MatrixType = _MatrixType
- using PermutationType = PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime>
- using PlainObject = MatrixType::PlainObject
- using TranspositionType = Transpositions<RowsAtCompileTime, MaxRowsAtCompileTime>
Constructors, destructors, conversion operators
- PartialPivLU()
- Default Constructor.
- PartialPivLU(Index size) explicit
- Default Constructor with memory preallocation.
-
template<typename InputType>PartialPivLU(const EigenBase<InputType>& matrix) explicit
-
template<typename InputType>PartialPivLU(EigenBase<InputType>& matrix) explicit
Public functions
- auto cols() const -> Index
-
template<typename InputType>auto compute(const EigenBase<InputType>& matrix) -> PartialPivLU&
- auto determinant() const -> Scalar
- auto inverse() const -> const Inverse<PartialPivLU>
- auto matrixLU() const -> const MatrixType&
- auto permutationP() const -> const PermutationType&
- auto rcond() const -> RealScalar
- auto reconstructedMatrix() const -> MatrixType
- auto rows() const -> Index
-
template<typename Rhs>auto solve(const MatrixBase<Rhs>& b) const -> const Solve<PartialPivLU, Rhs>
Protected static functions
- static void check_template_parameters()
Protected functions
- void compute()
Protected variables
- signed char m_det_p
- bool m_isInitialized
- RealScalar m_l1_norm
- MatrixType m_lu
- PermutationType m_p
- TranspositionType m_rowsTranspositions
Enum documentation
template<typename _MatrixType>
enum Eigen:: PartialPivLU<_MatrixType>:: (anonymous)
Enumerators | |
---|---|
MaxRowsAtCompileTime | |
MaxColsAtCompileTime |
Typedef documentation
template<typename _MatrixType>
typedef SolverBase<PartialPivLU> Eigen:: PartialPivLU<_MatrixType>:: Base
template<typename _MatrixType>
typedef _MatrixType Eigen:: PartialPivLU<_MatrixType>:: MatrixType
template<typename _MatrixType>
typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime> Eigen:: PartialPivLU<_MatrixType>:: PermutationType
template<typename _MatrixType>
typedef MatrixType::PlainObject Eigen:: PartialPivLU<_MatrixType>:: PlainObject
template<typename _MatrixType>
typedef Transpositions<RowsAtCompileTime, MaxRowsAtCompileTime> Eigen:: PartialPivLU<_MatrixType>:: TranspositionType
Function documentation
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&).
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.
template<typename _MatrixType>
template<typename InputType>
Eigen:: PartialPivLU<_MatrixType>:: PartialPivLU(const EigenBase<InputType>& matrix) explicit
Parameters | |
---|---|
matrix | the matrix of which to compute the LU decomposition. |
Constructor.
template<typename _MatrixType>
template<typename InputType>
Eigen:: PartialPivLU<_MatrixType>:: PartialPivLU(EigenBase<InputType>& matrix) explicit
Parameters | |
---|---|
matrix | the matrix of which to compute the LU decomposition. |
Constructor for inplace decomposition
template<typename _MatrixType>
template<typename InputType>
PartialPivLU& Eigen:: PartialPivLU<_MatrixType>:: compute(const EigenBase<InputType>& matrix)
template<typename _MatrixType>
Scalar Eigen:: PartialPivLU<_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. |
---|
template<typename _MatrixType>
const Inverse<PartialPivLU> Eigen:: PartialPivLU<_MatrixType>:: inverse() const
Returns | the inverse of the matrix of which *this is the LU decomposition. |
---|
template<typename _MatrixType>
const MatrixType& Eigen:: PartialPivLU<_MatrixType>:: matrixLU() const
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). |
---|
template<typename _MatrixType>
const PermutationType& Eigen:: PartialPivLU<_MatrixType>:: permutationP() const
Returns | the permutation matrix P. |
---|
template<typename _MatrixType>
RealScalar Eigen:: PartialPivLU<_MatrixType>:: rcond() const
Returns | an estimate of the reciprocal condition number of the matrix of which *this is the LU decomposition. |
---|
template<typename _MatrixType>
MatrixType Eigen:: PartialPivLU<_MatrixType>:: reconstructedMatrix() const
Returns | the matrix represented by the decomposition, i.e., it returns the product: P^{-1} L U. This function is provided for debug purpose. |
---|
template<typename _MatrixType>
template<typename Rhs>
const Solve<PartialPivLU, Rhs> Eigen:: PartialPivLU<_MatrixType>:: solve(const MatrixBase<Rhs>& b) const
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. |
Returns | the solution. |
This method returns the solution x to the equation Ax=b, where A is the matrix of which *this is the LU decomposition.
Example:
MatrixXd A = MatrixXd::Random(3,3); MatrixXd B = MatrixXd::Random(3,2); 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;
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.