![]() |
Eigen
3.4.90 (git rev 67eeba6e720c5745abc77ae6c92ce0a44aa7b7ae)
|
LU decomposition of a matrix with partial pivoting, and related features.
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.
Public Member Functions | |
Scalar | determinant () const |
const Inverse< PartialPivLU > | inverse () 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 PermutationType & | permutationP () const |
RealScalar | rcond () const |
MatrixType | reconstructedMatrix () const |
template<typename Rhs > | |
const Solve< PartialPivLU, Rhs > | solve (const MatrixBase< Rhs > &b) const |
![]() | |
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 |
![]() | |
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 | |
![]() | |
typedef Eigen::Index | Index |
The interface type of indices. More... | |
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&).
|
explicit |
Default Constructor with memory preallocation.
Like the default constructor but with preallocation of the internal data according to the specified problem size.
|
explicit |
Constructor.
matrix | the matrix of which to compute the LU decomposition. |
|
explicit |
Constructor for inplace decomposition
matrix | the matrix of which to compute the LU decomposition. |
PartialPivLU< MatrixType >::Scalar Eigen::PartialPivLU< MatrixType >::determinant |
|
inline |
|
inline |
|
inline |
|
inline |
*this
is the LU decomposition. MatrixType Eigen::PartialPivLU< MatrixType >::reconstructedMatrix |
|
inline |
This method returns the solution x to the equation Ax=b, where A is the matrix of which *this is the LU decomposition.
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. |
Example:
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.