Eigen::PartialPivLU< MatrixType_, PermutationIndex_ > Class Template Reference

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

+ Inheritance diagram for Eigen::PartialPivLU< MatrixType_, PermutationIndex_ >:

Public Types

enum  {
  MaxRowsAtCompileTime ,
  MaxColsAtCompileTime
}
 
typedef SolverBase< PartialPivLUBase
 
typedef MatrixType_ MatrixType
 
using PermutationIndex = PermutationIndex_
 
typedef PermutationMatrix< RowsAtCompileTime, MaxRowsAtCompileTime, PermutationIndexPermutationType
 
typedef MatrixType::PlainObject PlainObject
 
typedef Transpositions< RowsAtCompileTime, MaxRowsAtCompileTime, PermutationIndexTranspositionType
 
- Public Types inherited from Eigen::SolverBase< PartialPivLU< MatrixType_, PermutationIndex_ > >
enum  
 
typedef std::conditional_t< NumTraits< Scalar >::IsComplex, CwiseUnaryOp< internal::scalar_conjugate_op< Scalar >, const ConstTransposeReturnType >, const ConstTransposeReturnTypeAdjointReturnType
 
typedef EigenBase< PartialPivLU< MatrixType_, PermutationIndex_ > > Base
 
typedef Scalar CoeffReturnType
 
typedef Transpose< const PartialPivLU< MatrixType_, PermutationIndex_ > > ConstTransposeReturnType
 
typedef internal::traits< PartialPivLU< MatrixType_, PermutationIndex_ > >::Scalar Scalar
 
- Public Types inherited from Eigen::EigenBase< Derived >
typedef Eigen::Index Index
 The interface type of indices. More...
 
typedef internal::traits< Derived >::StorageKind StorageKind
 

Public Member Functions

EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
template<typename InputType >
PartialPivLUcompute (const EigenBase< InputType > &matrix)
 
Scalar determinant () const
 
const Inverse< PartialPivLUinverse () const
 
const MatrixTypematrixLU () 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
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
template<typename Rhs >
const Solve< PartialPivLU, Rhs > solve (const MatrixBase< Rhs > &b) const
 
- Public Member Functions inherited from Eigen::SolverBase< PartialPivLU< MatrixType_, PermutationIndex_ > >
const AdjointReturnType adjoint () const
 
PartialPivLU< MatrixType_, PermutationIndex_ > & derived ()
 
const PartialPivLU< MatrixType_, PermutationIndex_ > & derived () const
 
const Solve< PartialPivLU< MatrixType_, PermutationIndex_ >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
 SolverBase ()
 
const ConstTransposeReturnType transpose () const
 
 ~SolverBase ()
 
- Public Member Functions inherited from Eigen::EigenBase< Derived >
template<typename Dest >
void addTo (Dest &dst) const
 
template<typename Dest >
void applyThisOnTheLeft (Dest &dst) const
 
template<typename Dest >
void applyThisOnTheRight (Dest &dst) const
 
EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
Derived & const_cast_derived () const
 
const Derived & const_derived () const
 
Derived & derived ()
 
const Derived & derived () const
 
template<typename Dest >
void evalTo (Dest &dst) const
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
EIGEN_CONSTEXPR Index size () const EIGEN_NOEXCEPT
 
template<typename Dest >
void subTo (Dest &dst) const
 

Protected Member Functions

void compute ()
 
- Protected Member Functions inherited from Eigen::SolverBase< PartialPivLU< MatrixType_, PermutationIndex_ > >
void _check_solve_assertion (const Rhs &b) const
 

Protected Attributes

signed char m_det_p
 
bool m_isInitialized
 
RealScalar m_l1_norm
 
MatrixType m_lu
 
PermutationType m_p
 
TranspositionType m_rowsTranspositions
 

Detailed Description

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

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

Definition at line 78 of file PartialPivLU.h.

Member Typedef Documentation

◆ Base

template<typename MatrixType_ , typename PermutationIndex_ >
typedef SolverBase<PartialPivLU> Eigen::PartialPivLU< MatrixType_, PermutationIndex_ >::Base

Definition at line 84 of file PartialPivLU.h.

◆ MatrixType

template<typename MatrixType_ , typename PermutationIndex_ >
typedef MatrixType_ Eigen::PartialPivLU< MatrixType_, PermutationIndex_ >::MatrixType

Definition at line 83 of file PartialPivLU.h.

◆ PermutationIndex

template<typename MatrixType_ , typename PermutationIndex_ >
using Eigen::PartialPivLU< MatrixType_, PermutationIndex_ >::PermutationIndex = PermutationIndex_

Definition at line 92 of file PartialPivLU.h.

◆ PermutationType

template<typename MatrixType_ , typename PermutationIndex_ >
typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime, PermutationIndex> Eigen::PartialPivLU< MatrixType_, PermutationIndex_ >::PermutationType

Definition at line 93 of file PartialPivLU.h.

◆ PlainObject

template<typename MatrixType_ , typename PermutationIndex_ >
typedef MatrixType::PlainObject Eigen::PartialPivLU< MatrixType_, PermutationIndex_ >::PlainObject

Definition at line 95 of file PartialPivLU.h.

◆ TranspositionType

template<typename MatrixType_ , typename PermutationIndex_ >
typedef Transpositions<RowsAtCompileTime, MaxRowsAtCompileTime, PermutationIndex> Eigen::PartialPivLU< MatrixType_, PermutationIndex_ >::TranspositionType

Definition at line 94 of file PartialPivLU.h.

Member Enumeration Documentation

◆ anonymous enum

template<typename MatrixType_ , typename PermutationIndex_ >
anonymous enum
Enumerator
MaxRowsAtCompileTime 
MaxColsAtCompileTime 

Definition at line 88 of file PartialPivLU.h.

88  {
89  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
90  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
91  };

Constructor & Destructor Documentation

◆ PartialPivLU() [1/4]

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

Default Constructor.

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

Definition at line 284 of file PartialPivLU.h.

285  : m_lu(),
286  m_p(),
288  m_l1_norm(0),
289  m_det_p(0),
290  m_isInitialized(false)
291 {
292 }
RealScalar m_l1_norm
Definition: PartialPivLU.h:278
TranspositionType m_rowsTranspositions
Definition: PartialPivLU.h:277
PermutationType m_p
Definition: PartialPivLU.h:276
signed char m_det_p
Definition: PartialPivLU.h:279

◆ PartialPivLU() [2/4]

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

Definition at line 295 of file PartialPivLU.h.

296  : m_lu(size, size),
297  m_p(size),
299  m_l1_norm(0),
300  m_det_p(0),
301  m_isInitialized(false)
302 {
303 }
EIGEN_CONSTEXPR Index size() const EIGEN_NOEXCEPT
Definition: EigenBase.h:69

◆ PartialPivLU() [3/4]

template<typename MatrixType , typename PermutationIndex >
template<typename InputType >
Eigen::PartialPivLU< MatrixType, PermutationIndex >::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.

Definition at line 307 of file PartialPivLU.h.

308  : m_lu(matrix.rows(),matrix.cols()),
309  m_p(matrix.rows()),
310  m_rowsTranspositions(matrix.rows()),
311  m_l1_norm(0),
312  m_det_p(0),
313  m_isInitialized(false)
314 {
315  compute(matrix.derived());
316 }

◆ PartialPivLU() [4/4]

template<typename MatrixType , typename PermutationIndex >
template<typename InputType >
Eigen::PartialPivLU< MatrixType, PermutationIndex >::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.

Definition at line 320 of file PartialPivLU.h.

321  : m_lu(matrix.derived()),
322  m_p(matrix.rows()),
323  m_rowsTranspositions(matrix.rows()),
324  m_l1_norm(0),
325  m_det_p(0),
326  m_isInitialized(false)
327 {
328  compute();
329 }

Member Function Documentation

◆ cols()

template<typename MatrixType_ , typename PermutationIndex_ >
EIGEN_CONSTEXPR Index Eigen::PartialPivLU< MatrixType_, PermutationIndex_ >::cols ( void  ) const
inline

Definition at line 223 of file PartialPivLU.h.

223 { return m_lu.cols(); }

◆ compute() [1/2]

template<typename MatrixType , typename PermutationIndex >
void Eigen::PartialPivLU< MatrixType, PermutationIndex >::compute
protected

Definition at line 525 of file PartialPivLU.h.

526 {
527  eigen_assert(m_lu.rows()<NumTraits<PermutationIndex>::highest());
528 
529  if(m_lu.cols()>0)
530  m_l1_norm = m_lu.cwiseAbs().colwise().sum().maxCoeff();
531  else
532  m_l1_norm = RealScalar(0);
533 
534  eigen_assert(m_lu.rows() == m_lu.cols() && "PartialPivLU is only for square (and moreover invertible) matrices");
535  const Index size = m_lu.rows();
536 
538 
539  typename TranspositionType::StorageIndex nb_transpositions;
541  m_det_p = (nb_transpositions%2) ? -1 : 1;
542 
544 
545  m_isInitialized = true;
546 }
#define eigen_assert(x)
Definition: Macros.h:902
void resize(Index newSize)
void partial_lu_inplace(MatrixType &lu, TranspositionType &row_transpositions, typename TranspositionType::StorageIndex &nb_transpositions)
Definition: PartialPivLU.h:505
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:41

◆ compute() [2/2]

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

Definition at line 134 of file PartialPivLU.h.

134  {
135  m_lu = matrix.derived();
136  compute();
137  return *this;
138  }

◆ determinant()

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

Definition at line 549 of file PartialPivLU.h.

550 {
551  eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
552  return Scalar(m_det_p) * m_lu.diagonal().prod();
553 }
internal::traits< PartialPivLU< MatrixType_, PermutationIndex_ > >::Scalar Scalar
Definition: SolverBase.h:75

◆ inverse()

template<typename MatrixType_ , typename PermutationIndex_ >
const Inverse<PartialPivLU> Eigen::PartialPivLU< MatrixType_, PermutationIndex_ >::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()

Definition at line 199 of file PartialPivLU.h.

200  {
201  eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
202  return Inverse<PartialPivLU>(*this);
203  }

◆ matrixLU()

template<typename MatrixType_ , typename PermutationIndex_ >
const MatrixType& Eigen::PartialPivLU< MatrixType_, PermutationIndex_ >::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 146 of file PartialPivLU.h.

147  {
148  eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
149  return m_lu;
150  }

◆ permutationP()

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

Definition at line 154 of file PartialPivLU.h.

155  {
156  eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
157  return m_p;
158  }

◆ rcond()

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

Definition at line 186 of file PartialPivLU.h.

187  {
188  eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
190  }
Decomposition::RealScalar rcond_estimate_helper(typename Decomposition::RealScalar matrix_norm, const Decomposition &dec)
Reciprocal condition number estimator.

◆ reconstructedMatrix()

template<typename MatrixType , typename PermutationIndex >
MatrixType Eigen::PartialPivLU< MatrixType, PermutationIndex >::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.

Definition at line 559 of file PartialPivLU.h.

560 {
561  eigen_assert(m_isInitialized && "LU is not initialized.");
562  // LU
563  MatrixType res = m_lu.template triangularView<UnitLower>().toDenseMatrix()
564  * m_lu.template triangularView<Upper>();
565 
566  // P^{-1}(LU)
567  res = m_p.inverse() * res;
568 
569  return res;
570 }
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
Matrix< float, 1, Dynamic > MatrixType
InverseReturnType inverse() const

◆ rows()

template<typename MatrixType_ , typename PermutationIndex_ >
EIGEN_CONSTEXPR Index Eigen::PartialPivLU< MatrixType_, PermutationIndex_ >::rows ( void  ) const
inline

Definition at line 222 of file PartialPivLU.h.

222 { return m_lu.rows(); }

◆ solve()

template<typename MatrixType_ , typename PermutationIndex_ >
template<typename Rhs >
const Solve<PartialPivLU, Rhs> Eigen::PartialPivLU< MatrixType_, PermutationIndex_ >::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;
MatrixXcf A
MatrixXf B
static const RandomReturnType Random()
Definition: Random.h:114
Matrix< double, Dynamic, Dynamic > MatrixXd
Dynamic×Dynamic matrix of type double.
Definition: Matrix.h:502

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()

Member Data Documentation

◆ m_det_p

template<typename MatrixType_ , typename PermutationIndex_ >
signed char Eigen::PartialPivLU< MatrixType_, PermutationIndex_ >::m_det_p
protected

Definition at line 279 of file PartialPivLU.h.

◆ m_isInitialized

template<typename MatrixType_ , typename PermutationIndex_ >
bool Eigen::PartialPivLU< MatrixType_, PermutationIndex_ >::m_isInitialized
protected

Definition at line 280 of file PartialPivLU.h.

◆ m_l1_norm

template<typename MatrixType_ , typename PermutationIndex_ >
RealScalar Eigen::PartialPivLU< MatrixType_, PermutationIndex_ >::m_l1_norm
protected

Definition at line 278 of file PartialPivLU.h.

◆ m_lu

template<typename MatrixType_ , typename PermutationIndex_ >
MatrixType Eigen::PartialPivLU< MatrixType_, PermutationIndex_ >::m_lu
protected

Definition at line 275 of file PartialPivLU.h.

◆ m_p

template<typename MatrixType_ , typename PermutationIndex_ >
PermutationType Eigen::PartialPivLU< MatrixType_, PermutationIndex_ >::m_p
protected

Definition at line 276 of file PartialPivLU.h.

◆ m_rowsTranspositions

template<typename MatrixType_ , typename PermutationIndex_ >
TranspositionType Eigen::PartialPivLU< MatrixType_, PermutationIndex_ >::m_rowsTranspositions
protected

Definition at line 277 of file PartialPivLU.h.


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