Public Types | Public Member Functions | Static Protected 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>

List of all members.

Public Types

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

Public Member Functions

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

Static Protected Member Functions

static void check_template_parameters ()

Protected Attributes

Index m_det_p
bool m_isInitialized
MatrixType m_lu
PermutationType m_p
TranspositionType m_rowsTranspositions

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


Member Typedef Documentation

template<typename _MatrixType>
typedef MatrixType::Index Eigen::PartialPivLU< _MatrixType >::Index

Definition at line 62 of file PartialPivLU.h.

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

Definition at line 51 of file PartialPivLU.h.

template<typename _MatrixType>
typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime> Eigen::PartialPivLU< _MatrixType >::PermutationType

Definition at line 63 of file PartialPivLU.h.

template<typename _MatrixType>
typedef NumTraits<typename MatrixType::Scalar>::Real Eigen::PartialPivLU< _MatrixType >::RealScalar

Definition at line 60 of file PartialPivLU.h.

template<typename _MatrixType>
typedef MatrixType::Scalar Eigen::PartialPivLU< _MatrixType >::Scalar

Definition at line 59 of file PartialPivLU.h.

template<typename _MatrixType>
typedef internal::traits<MatrixType>::StorageKind Eigen::PartialPivLU< _MatrixType >::StorageKind

Definition at line 61 of file PartialPivLU.h.

template<typename _MatrixType>
typedef Transpositions<RowsAtCompileTime, MaxRowsAtCompileTime> Eigen::PartialPivLU< _MatrixType >::TranspositionType

Definition at line 64 of file PartialPivLU.h.


Member Enumeration Documentation

template<typename _MatrixType>
anonymous enum
Enumerator:
RowsAtCompileTime 
ColsAtCompileTime 
Options 
MaxRowsAtCompileTime 
MaxColsAtCompileTime 

Definition at line 52 of file PartialPivLU.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 188 of file PartialPivLU.h.

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

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


Member Function Documentation

template<typename _MatrixType>
static void Eigen::PartialPivLU< _MatrixType >::check_template_parameters ( ) [inline, static, protected]

Definition at line 175 of file PartialPivLU.h.

template<typename _MatrixType>
Index Eigen::PartialPivLU< _MatrixType >::cols ( void  ) const [inline]

Definition at line 171 of file PartialPivLU.h.

template<typename MatrixType >
PartialPivLU< MatrixType > & Eigen::PartialPivLU< MatrixType >::compute ( const MatrixType matrix)

Definition at line 393 of file PartialPivLU.h.

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

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.

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.

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

Definition at line 108 of file PartialPivLU.h.

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

template<typename _MatrixType>
Index Eigen::PartialPivLU< _MatrixType >::rows ( void  ) const [inline]

Definition at line 170 of file PartialPivLU.h.

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.


Member Data Documentation

template<typename _MatrixType>
Index Eigen::PartialPivLU< _MatrixType >::m_det_p [protected]

Definition at line 183 of file PartialPivLU.h.

template<typename _MatrixType>
bool Eigen::PartialPivLU< _MatrixType >::m_isInitialized [protected]

Definition at line 184 of file PartialPivLU.h.

template<typename _MatrixType>
MatrixType Eigen::PartialPivLU< _MatrixType >::m_lu [protected]

Definition at line 180 of file PartialPivLU.h.

template<typename _MatrixType>
PermutationType Eigen::PartialPivLU< _MatrixType >::m_p [protected]

Definition at line 181 of file PartialPivLU.h.

template<typename _MatrixType>
TranspositionType Eigen::PartialPivLU< _MatrixType >::m_rowsTranspositions [protected]

Definition at line 182 of file PartialPivLU.h.


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


turtlebot_exploration_3d
Author(s): Bona , Shawn
autogenerated on Thu Jun 6 2019 21:00:53