Public Types | Public Member Functions | Protected Member Functions | Static Protected Member Functions | Protected Attributes | List of all members
Eigen::FullPivLU< _MatrixType > Class Template Reference

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

#include <ForwardDeclarations.h>

Inheritance diagram for Eigen::FullPivLU< _MatrixType >:
Inheritance graph
[legend]

Public Types

enum  { MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime }
 
typedef SolverBase< FullPivLUBase
 
typedef internal::plain_col_type< MatrixType, StorageIndex >::type IntColVectorType
 
typedef internal::plain_row_type< MatrixType, StorageIndex >::type IntRowVectorType
 
typedef _MatrixType MatrixType
 
typedef PermutationMatrix< RowsAtCompileTime, MaxRowsAtCompileTimePermutationPType
 
typedef PermutationMatrix< ColsAtCompileTime, MaxColsAtCompileTimePermutationQType
 
typedef MatrixType::PlainObject PlainObject
 
- Public Types inherited from Eigen::SolverBase< FullPivLU< _MatrixType > >
enum  
 
typedef internal::conditional< NumTraits< Scalar >::IsComplex, CwiseUnaryOp< internal::scalar_conjugate_op< Scalar >, ConstTransposeReturnType >, ConstTransposeReturnType >::type AdjointReturnType
 
typedef EigenBase< FullPivLU< _MatrixType > > Base
 
typedef Scalar CoeffReturnType
 
typedef internal::add_const< Transpose< const FullPivLU< _MatrixType > > >::type ConstTransposeReturnType
 
typedef internal::traits< FullPivLU< _MatrixType > >::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

template<typename RhsType , typename DstType >
EIGEN_DEVICE_FUNC void _solve_impl (const RhsType &rhs, DstType &dst) const
 
template<typename RhsType , typename DstType >
void _solve_impl (const RhsType &rhs, DstType &dst) const
 
template<bool Conjugate, typename RhsType , typename DstType >
EIGEN_DEVICE_FUNC void _solve_impl_transposed (const RhsType &rhs, DstType &dst) const
 
template<bool Conjugate, typename RhsType , typename DstType >
void _solve_impl_transposed (const RhsType &rhs, DstType &dst) const
 
EIGEN_DEVICE_FUNC Index cols () const
 
template<typename InputType >
FullPivLUcompute (const EigenBase< InputType > &matrix)
 
internal::traits< MatrixType >::Scalar determinant () const
 
Index dimensionOfKernel () const
 
 FullPivLU ()
 Default Constructor. More...
 
 FullPivLU (Index rows, Index cols)
 Default Constructor with memory preallocation. More...
 
template<typename InputType >
 FullPivLU (const EigenBase< InputType > &matrix)
 
template<typename InputType >
 FullPivLU (EigenBase< InputType > &matrix)
 Constructs a LU factorization from a given matrix. More...
 
const internal::image_retval< FullPivLUimage (const MatrixType &originalMatrix) const
 
const Inverse< FullPivLUinverse () const
 
bool isInjective () const
 
bool isInvertible () const
 
bool isSurjective () const
 
const internal::kernel_retval< FullPivLUkernel () const
 
const MatrixTypematrixLU () const
 
RealScalar maxPivot () const
 
Index nonzeroPivots () const
 
EIGEN_DEVICE_FUNC const PermutationPTypepermutationP () const
 
const PermutationQTypepermutationQ () const
 
Index rank () const
 
RealScalar rcond () const
 
MatrixType reconstructedMatrix () const
 
EIGEN_DEVICE_FUNC Index rows () const
 
FullPivLUsetThreshold (const RealScalar &threshold)
 
FullPivLUsetThreshold (Default_t)
 
template<typename Rhs >
const Solve< FullPivLU, Rhs > solve (const MatrixBase< Rhs > &b) const
 
RealScalar threshold () const
 
- Public Member Functions inherited from Eigen::SolverBase< FullPivLU< _MatrixType > >
AdjointReturnType adjoint () const
 
const Solve< FullPivLU< _MatrixType >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
 SolverBase ()
 
ConstTransposeReturnType transpose () const
 
 ~SolverBase ()
 
- Public Member Functions inherited from Eigen::EigenBase< Derived >
template<typename Dest >
EIGEN_DEVICE_FUNC void addTo (Dest &dst) const
 
template<typename Dest >
EIGEN_DEVICE_FUNC void applyThisOnTheLeft (Dest &dst) const
 
template<typename Dest >
EIGEN_DEVICE_FUNC void applyThisOnTheRight (Dest &dst) const
 
EIGEN_DEVICE_FUNC Index cols () const
 
EIGEN_DEVICE_FUNC Derived & const_cast_derived () const
 
EIGEN_DEVICE_FUNC const Derived & const_derived () const
 
EIGEN_DEVICE_FUNC Derived & derived ()
 
EIGEN_DEVICE_FUNC const Derived & derived () const
 
template<typename Dest >
EIGEN_DEVICE_FUNC void evalTo (Dest &dst) const
 
EIGEN_DEVICE_FUNC Index rows () const
 
EIGEN_DEVICE_FUNC Index size () const
 
template<typename Dest >
EIGEN_DEVICE_FUNC void subTo (Dest &dst) const
 

Protected Member Functions

void computeInPlace ()
 

Static Protected Member Functions

static void check_template_parameters ()
 

Protected Attributes

IntRowVectorType m_colsTranspositions
 
signed char m_det_pq
 
bool m_isInitialized
 
RealScalar m_l1_norm
 
MatrixType m_lu
 
RealScalar m_maxpivot
 
Index m_nonzero_pivots
 
PermutationPType m_p
 
RealScalar m_prescribedThreshold
 
PermutationQType m_q
 
IntColVectorType m_rowsTranspositions
 
bool m_usePrescribedThreshold
 

Detailed Description

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

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

Template Parameters
_MatrixTypethe type of the matrix of which we are computing the LU decomposition

This class represents a LU decomposition of any matrix, with complete pivoting: the matrix A is decomposed as $ A = P^{-1} L U Q^{-1} $ where L is unit-lower-triangular, U is upper-triangular, and P and Q are permutation matrices. This is a rank-revealing LU decomposition. The eigenvalues (diagonal coefficients) of U are sorted in such a way that any zeros are at the end.

This decomposition provides the generic approach to solving systems of linear equations, computing the rank, invertibility, inverse, kernel, and determinant.

This LU decomposition is very stable and well tested with large matrices. However there are use cases where the SVD decomposition is inherently more stable and/or flexible. For example, when computing the kernel of a matrix, working with the SVD allows to select the smallest singular values of the matrix, something that the LU decomposition doesn't see.

The data of the LU decomposition can be directly accessed through the methods matrixLU(), permutationP(), permutationQ().

As an exemple, here is how the original matrix can be retrieved:

typedef Matrix<double, 5, 3> Matrix5x3;
typedef Matrix<double, 5, 5> Matrix5x5;
Matrix5x3 m = Matrix5x3::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is, up to permutations, its LU decomposition matrix:"
<< endl << lu.matrixLU() << endl;
cout << "Here is the L part:" << endl;
Matrix5x5 l = Matrix5x5::Identity();
l.block<5,3>(0,0).triangularView<StrictlyLower>() = lu.matrixLU();
cout << l << endl;
cout << "Here is the U part:" << endl;
Matrix5x3 u = lu.matrixLU().triangularView<Upper>();
cout << u << endl;
cout << "Let us now reconstruct the original matrix m:" << endl;
cout << lu.permutationP().inverse() * l * u * lu.permutationQ().inverse() << endl;

Output:

This class supports the inplace decomposition mechanism.

See also
MatrixBase::fullPivLu(), MatrixBase::determinant(), MatrixBase::inverse()

Definition at line 249 of file ForwardDeclarations.h.

Member Typedef Documentation

template<typename _MatrixType>
typedef SolverBase<FullPivLU> Eigen::FullPivLU< _MatrixType >::Base

Definition at line 64 of file FullPivLU.h.

template<typename _MatrixType>
typedef internal::plain_col_type<MatrixType, StorageIndex>::type Eigen::FullPivLU< _MatrixType >::IntColVectorType

Definition at line 73 of file FullPivLU.h.

template<typename _MatrixType>
typedef internal::plain_row_type<MatrixType, StorageIndex>::type Eigen::FullPivLU< _MatrixType >::IntRowVectorType

Definition at line 72 of file FullPivLU.h.

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

Definition at line 63 of file FullPivLU.h.

template<typename _MatrixType>
typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime> Eigen::FullPivLU< _MatrixType >::PermutationPType

Definition at line 75 of file FullPivLU.h.

template<typename _MatrixType>
typedef PermutationMatrix<ColsAtCompileTime, MaxColsAtCompileTime> Eigen::FullPivLU< _MatrixType >::PermutationQType

Definition at line 74 of file FullPivLU.h.

template<typename _MatrixType>
typedef MatrixType::PlainObject Eigen::FullPivLU< _MatrixType >::PlainObject

Definition at line 76 of file FullPivLU.h.

Member Enumeration Documentation

template<typename _MatrixType>
anonymous enum
Enumerator
MaxRowsAtCompileTime 
MaxColsAtCompileTime 

Definition at line 68 of file FullPivLU.h.

Constructor & Destructor Documentation

template<typename MatrixType >
Eigen::FullPivLU< MatrixType >::FullPivLU ( )

Default Constructor.

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

Definition at line 444 of file FullPivLU.h.

template<typename MatrixType >
Eigen::FullPivLU< MatrixType >::FullPivLU ( Index  rows,
Index  cols 
)

Default Constructor with memory preallocation.

Like the default constructor but with preallocation of the internal data according to the specified problem size.

See also
FullPivLU()

Definition at line 450 of file FullPivLU.h.

template<typename MatrixType >
template<typename InputType >
Eigen::FullPivLU< MatrixType >::FullPivLU ( const EigenBase< InputType > &  matrix)
explicit

Constructor.

Parameters
matrixthe matrix of which to compute the LU decomposition. It is required to be nonzero.

Definition at line 463 of file FullPivLU.h.

template<typename MatrixType >
template<typename InputType >
Eigen::FullPivLU< MatrixType >::FullPivLU ( EigenBase< InputType > &  matrix)
explicit

Constructs a LU factorization from a given matrix.

This overloaded constructor is provided for inplace decomposition when MatrixType is a Eigen::Ref.

See also
FullPivLU(const EigenBase&)

Definition at line 477 of file FullPivLU.h.

Member Function Documentation

template<typename _MatrixType>
template<typename RhsType , typename DstType >
EIGEN_DEVICE_FUNC void Eigen::FullPivLU< _MatrixType >::_solve_impl ( const RhsType &  rhs,
DstType &  dst 
) const
template<typename _MatrixType>
template<typename RhsType , typename DstType >
void Eigen::FullPivLU< _MatrixType >::_solve_impl ( const RhsType &  rhs,
DstType &  dst 
) const

Definition at line 747 of file FullPivLU.h.

template<typename _MatrixType>
template<bool Conjugate, typename RhsType , typename DstType >
EIGEN_DEVICE_FUNC void Eigen::FullPivLU< _MatrixType >::_solve_impl_transposed ( const RhsType &  rhs,
DstType &  dst 
) const
template<typename _MatrixType>
template<bool Conjugate, typename RhsType , typename DstType >
void Eigen::FullPivLU< _MatrixType >::_solve_impl_transposed ( const RhsType &  rhs,
DstType &  dst 
) const

Definition at line 795 of file FullPivLU.h.

template<typename _MatrixType>
static void Eigen::FullPivLU< _MatrixType >::check_template_parameters ( )
inlinestaticprotected

Definition at line 424 of file FullPivLU.h.

template<typename _MatrixType>
EIGEN_DEVICE_FUNC Index Eigen::FullPivLU< _MatrixType >::cols ( void  ) const
inline

Definition at line 410 of file FullPivLU.h.

template<typename _MatrixType>
template<typename InputType >
FullPivLU& Eigen::FullPivLU< _MatrixType >::compute ( const EigenBase< InputType > &  matrix)
inline

Computes the LU decomposition of the given matrix.

Parameters
matrixthe matrix of which to compute the LU decomposition. It is required to be nonzero.
Returns
a reference to *this

Definition at line 119 of file FullPivLU.h.

template<typename MatrixType >
void Eigen::FullPivLU< MatrixType >::computeInPlace ( )
protected

Definition at line 490 of file FullPivLU.h.

template<typename MatrixType >
internal::traits< MatrixType >::Scalar Eigen::FullPivLU< 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
This is only for square matrices.
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 583 of file FullPivLU.h.

template<typename _MatrixType>
Index Eigen::FullPivLU< _MatrixType >::dimensionOfKernel ( ) const
inline
Returns
the dimension of the kernel of the matrix of which *this is the LU decomposition.
Note
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).

Definition at line 349 of file FullPivLU.h.

template<typename _MatrixType>
const internal::image_retval<FullPivLU> Eigen::FullPivLU< _MatrixType >::image ( const MatrixType originalMatrix) const
inline
Returns
the image of the matrix, also called its column-space. The columns of the returned matrix will form a basis of the image (column-space).
Parameters
originalMatrixthe original matrix, of which *this is the LU decomposition. The reason why it is needed to pass it here, is that this allows a large optimization, as otherwise this method would need to reconstruct it from the LU decomposition.
Note
If the image has dimension zero, then the returned matrix is a column-vector filled with zeros.
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).

Example:

Matrix3d m;
m << 1,1,0,
1,3,2,
0,1,1;
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Notice that the middle column is the sum of the two others, so the "
<< "columns are linearly dependent." << endl;
cout << "Here is a matrix whose columns have the same span but are linearly independent:"
<< endl << m.fullPivLu().image(m) << endl;

Output:

See also
kernel()

Definition at line 215 of file FullPivLU.h.

template<typename _MatrixType>
const Inverse<FullPivLU> Eigen::FullPivLU< _MatrixType >::inverse ( ) const
inline
Returns
the inverse of the matrix of which *this is the LU decomposition.
Note
If this matrix is not invertible, the returned matrix has undefined coefficients. Use isInvertible() to first determine whether this matrix is invertible.
See also
MatrixBase::inverse()

Definition at line 400 of file FullPivLU.h.

template<typename _MatrixType>
bool Eigen::FullPivLU< _MatrixType >::isInjective ( ) const
inline
Returns
true if the matrix of which *this is the LU decomposition represents an injective linear map, i.e. has trivial kernel; false otherwise.
Note
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).

Definition at line 362 of file FullPivLU.h.

template<typename _MatrixType>
bool Eigen::FullPivLU< _MatrixType >::isInvertible ( ) const
inline
Returns
true if the matrix of which *this is the LU decomposition is invertible.
Note
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).

Definition at line 387 of file FullPivLU.h.

template<typename _MatrixType>
bool Eigen::FullPivLU< _MatrixType >::isSurjective ( ) const
inline
Returns
true if the matrix of which *this is the LU decomposition represents a surjective linear map; false otherwise.
Note
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).

Definition at line 375 of file FullPivLU.h.

template<typename _MatrixType>
const internal::kernel_retval<FullPivLU> Eigen::FullPivLU< _MatrixType >::kernel ( ) const
inline
Returns
the kernel of the matrix, also called its null-space. The columns of the returned matrix will form a basis of the kernel.
Note
If the kernel has dimension zero, then the returned matrix is a column-vector filled with zeros.
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).

Example:

MatrixXf m = MatrixXf::Random(3,5);
cout << "Here is the matrix m:" << endl << m << endl;
MatrixXf ker = m.fullPivLu().kernel();
cout << "Here is a matrix whose columns form a basis of the kernel of m:"
<< endl << ker << endl;
cout << "By definition of the kernel, m*ker is zero:"
<< endl << m*ker << endl;

Output:

See also
image()

Definition at line 189 of file FullPivLU.h.

template<typename _MatrixType>
const MatrixType& Eigen::FullPivLU< _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 131 of file FullPivLU.h.

template<typename _MatrixType>
RealScalar Eigen::FullPivLU< _MatrixType >::maxPivot ( ) const
inline
Returns
the absolute value of the biggest pivot, i.e. the biggest diagonal coefficient of U.

Definition at line 153 of file FullPivLU.h.

template<typename _MatrixType>
Index Eigen::FullPivLU< _MatrixType >::nonzeroPivots ( ) const
inline
Returns
the number of nonzero pivots in the LU decomposition. Here nonzero is meant in the exact sense, not in a fuzzy sense. So that notion isn't really intrinsically interesting, but it is still useful when implementing algorithms.
See also
rank()

Definition at line 144 of file FullPivLU.h.

template<typename _MatrixType>
EIGEN_DEVICE_FUNC const PermutationPType& Eigen::FullPivLU< _MatrixType >::permutationP ( ) const
inline
Returns
the permutation matrix P
See also
permutationQ()

Definition at line 159 of file FullPivLU.h.

template<typename _MatrixType>
const PermutationQType& Eigen::FullPivLU< _MatrixType >::permutationQ ( ) const
inline
Returns
the permutation matrix Q
See also
permutationP()

Definition at line 169 of file FullPivLU.h.

template<typename _MatrixType>
Index Eigen::FullPivLU< _MatrixType >::rank ( ) const
inline
Returns
the rank of the matrix of which *this is the LU decomposition.
Note
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).

Definition at line 332 of file FullPivLU.h.

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

Definition at line 252 of file FullPivLU.h.

template<typename MatrixType >
MatrixType Eigen::FullPivLU< MatrixType >::reconstructedMatrix ( ) const
Returns
the matrix represented by the decomposition, i.e., it returns the product: $ P^{-1} L U Q^{-1} $. This function is provided for debug purposes.

Definition at line 594 of file FullPivLU.h.

template<typename _MatrixType>
EIGEN_DEVICE_FUNC Index Eigen::FullPivLU< _MatrixType >::rows ( void  ) const
inline

Definition at line 409 of file FullPivLU.h.

template<typename _MatrixType>
FullPivLU& Eigen::FullPivLU< _MatrixType >::setThreshold ( const RealScalar threshold)
inline

Allows to prescribe a threshold to be used by certain methods, such as rank(), who need to determine when pivots are to be considered nonzero. This is not used for the LU decomposition itself.

When it needs to get the threshold value, Eigen calls threshold(). By default, this uses a formula to automatically determine a reasonable threshold. Once you have called the present method setThreshold(const RealScalar&), your value is used instead.

Parameters
thresholdThe new value to use as the threshold.

A pivot will be considered nonzero if its absolute value is strictly greater than $ \vert pivot \vert \leqslant threshold \times \vert maxpivot \vert $ where maxpivot is the biggest pivot.

If you want to come back to the default behavior, call setThreshold(Default_t)

Definition at line 292 of file FullPivLU.h.

template<typename _MatrixType>
FullPivLU& Eigen::FullPivLU< _MatrixType >::setThreshold ( Default_t  )
inline

Allows to come back to the default behavior, letting Eigen use its default formula for determining the threshold.

You should pass the special object Eigen::Default as parameter here.

lu.setThreshold(Eigen::Default);

See the documentation of setThreshold(const RealScalar&).

Definition at line 307 of file FullPivLU.h.

template<typename _MatrixType>
template<typename Rhs >
const Solve<FullPivLU, Rhs> Eigen::FullPivLU< _MatrixType >::solve ( const MatrixBase< Rhs > &  b) const
inline
Returns
a 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
a solution.

Example:

Matrix<float,2,3> m = Matrix<float,2,3>::Random();
Matrix2f y = Matrix2f::Random();
cout << "Here is the matrix m:" << endl << m << endl;
cout << "Here is the matrix y:" << endl << y << endl;
Matrix<float,3,2> x = m.fullPivLu().solve(y);
if((m*x).isApprox(y))
{
cout << "Here is a solution x to the equation mx=y:" << endl << x << endl;
}
else
cout << "The equation mx=y does not have any solution." << endl;

Output:

See also
TriangularView::solve(), kernel(), inverse()

Definition at line 243 of file FullPivLU.h.

template<typename _MatrixType>
RealScalar Eigen::FullPivLU< _MatrixType >::threshold ( ) const
inline

Returns the threshold that will be used by certain methods such as rank().

See the documentation of setThreshold(const RealScalar&).

Definition at line 317 of file FullPivLU.h.

Member Data Documentation

template<typename _MatrixType>
IntRowVectorType Eigen::FullPivLU< _MatrixType >::m_colsTranspositions
protected

Definition at line 435 of file FullPivLU.h.

template<typename _MatrixType>
signed char Eigen::FullPivLU< _MatrixType >::m_det_pq
protected

Definition at line 439 of file FullPivLU.h.

template<typename _MatrixType>
bool Eigen::FullPivLU< _MatrixType >::m_isInitialized
protected

Definition at line 440 of file FullPivLU.h.

template<typename _MatrixType>
RealScalar Eigen::FullPivLU< _MatrixType >::m_l1_norm
protected

Definition at line 437 of file FullPivLU.h.

template<typename _MatrixType>
MatrixType Eigen::FullPivLU< _MatrixType >::m_lu
protected

Definition at line 431 of file FullPivLU.h.

template<typename _MatrixType>
RealScalar Eigen::FullPivLU< _MatrixType >::m_maxpivot
protected

Definition at line 438 of file FullPivLU.h.

template<typename _MatrixType>
Index Eigen::FullPivLU< _MatrixType >::m_nonzero_pivots
protected

Definition at line 436 of file FullPivLU.h.

template<typename _MatrixType>
PermutationPType Eigen::FullPivLU< _MatrixType >::m_p
protected

Definition at line 432 of file FullPivLU.h.

template<typename _MatrixType>
RealScalar Eigen::FullPivLU< _MatrixType >::m_prescribedThreshold
protected

Definition at line 438 of file FullPivLU.h.

template<typename _MatrixType>
PermutationQType Eigen::FullPivLU< _MatrixType >::m_q
protected

Definition at line 433 of file FullPivLU.h.

template<typename _MatrixType>
IntColVectorType Eigen::FullPivLU< _MatrixType >::m_rowsTranspositions
protected

Definition at line 434 of file FullPivLU.h.

template<typename _MatrixType>
bool Eigen::FullPivLU< _MatrixType >::m_usePrescribedThreshold
protected

Definition at line 440 of file FullPivLU.h.


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


gtsam
Author(s):
autogenerated on Sat May 8 2021 02:52:36