Shapeworks Studio  2.1
Shape analysis software suite
List of all members | Public Types | Public Member Functions | Protected Attributes
Eigen::PartialPivLU< _MatrixType > Class Template Reference

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

#include <PartialPivLU.h>

+ Collaboration diagram for Eigen::PartialPivLU< _MatrixType >:

Public Types

enum  {
  RowsAtCompileTime = MatrixType::RowsAtCompileTime, ColsAtCompileTime = MatrixType::ColsAtCompileTime, Options = MatrixType::Options, MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
}
 
typedef _MatrixType MatrixType
 
typedef MatrixType::Scalar Scalar
 
typedef NumTraits< typename MatrixType::Scalar >::Real RealScalar
 
typedef internal::traits< MatrixType >::StorageKind StorageKind
 
typedef MatrixType::Index Index
 
typedef PermutationMatrix< RowsAtCompileTime, MaxRowsAtCompileTime > PermutationType
 
typedef Transpositions< RowsAtCompileTime, MaxRowsAtCompileTime > TranspositionType
 

Public Member Functions

 PartialPivLU ()
 Default Constructor. More...
 
 PartialPivLU (Index size)
 Default Constructor with memory preallocation. More...
 
 PartialPivLU (const MatrixType &matrix)
 
PartialPivLUcompute (const MatrixType &matrix)
 
const MatrixType & matrixLU () const
 
const PermutationTypepermutationP () const
 
template<typename Rhs >
const internal::solve_retval< PartialPivLU, Rhs > solve (const MatrixBase< Rhs > &b) const
 
const internal::solve_retval< PartialPivLU, typename MatrixType::IdentityReturnType > inverse () const
 
internal::traits< MatrixType >::Scalar determinant () const
 
MatrixType reconstructedMatrix () const
 
Index rows () const
 
Index cols () const
 

Protected Attributes

MatrixType m_lu
 
PermutationType m_p
 
TranspositionType m_rowsTranspositions
 
Index m_det_p
 
bool m_isInitialized
 

Detailed Description

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

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

Parameters
MatrixTypethe 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().

See also
MatrixBase::partialPivLu(), MatrixBase::determinant(), MatrixBase::inverse(), MatrixBase::computeInverse(), class FullPivLU

Definition at line 217 of file ForwardDeclarations.h.

Constructor & Destructor 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&).

Definition at line 182 of file PartialPivLU.h.

183  : m_lu(),
184  m_p(),
185  m_rowsTranspositions(),
186  m_det_p(0),
187  m_isInitialized(false)
188 {
189 }
template<typename MatrixType >
Eigen::PartialPivLU< MatrixType >::PartialPivLU ( Index  size)

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 192 of file PartialPivLU.h.

193  : m_lu(size, size),
194  m_p(size),
195  m_rowsTranspositions(size),
196  m_det_p(0),
197  m_isInitialized(false)
198 {
199 }
template<typename MatrixType >
Eigen::PartialPivLU< MatrixType >::PartialPivLU ( const MatrixType &  matrix)

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 202 of file PartialPivLU.h.

203  : m_lu(matrix.rows(), matrix.rows()),
204  m_p(matrix.rows()),
205  m_rowsTranspositions(matrix.rows()),
206  m_det_p(0),
207  m_isInitialized(false)
208 {
209  compute(matrix);
210 }
Definition: math3d.h:219

Member Function Documentation

template<typename MatrixType >
internal::traits< 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.
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 410 of file PartialPivLU.h.

411 {
412  eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
413  return Scalar(m_det_p) * m_lu.diagonal().prod();
414 }
template<typename _MatrixType>
const internal::solve_retval<PartialPivLU,typename MatrixType::IdentityReturnType> Eigen::PartialPivLU< _MatrixType >::inverse ( void  ) 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 146 of file PartialPivLU.h.

147  {
148  eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
149  return internal::solve_retval<PartialPivLU,typename MatrixType::IdentityReturnType>
150  (*this, MatrixType::Identity(m_lu.rows(), m_lu.cols()));
151  }
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()

Definition at line 100 of file PartialPivLU.h.

101  {
102  eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
103  return m_lu;
104  }
template<typename _MatrixType>
const PermutationType& Eigen::PartialPivLU< _MatrixType >::permutationP ( ) const
inline
Returns
the permutation matrix P.

Definition at line 108 of file PartialPivLU.h.

109  {
110  eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
111  return m_p;
112  }
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.

Definition at line 420 of file PartialPivLU.h.

421 {
422  eigen_assert(m_isInitialized && "LU is not initialized.");
423  // LU
424  MatrixType res = m_lu.template triangularView<UnitLower>().toDenseMatrix()
425  * m_lu.template triangularView<Upper>();
426 
427  // P^{-1}(LU)
428  res = m_p.inverse() * res;
429 
430  return res;
431 }
Transpose< PermutationBase > inverse() const
template<typename _MatrixType>
template<typename Rhs >
const internal::solve_retval<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:

Output:

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

Definition at line 133 of file PartialPivLU.h.

134  {
135  eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
136  return internal::solve_retval<PartialPivLU, Rhs>(*this, b.derived());
137  }

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