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

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

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

Public Types

enum  {
  MaxRowsAtCompileTime ,
  MaxColsAtCompileTime
}
 
typedef SolverBase< FullPivLUBase
 
typedef internal::plain_col_type< MatrixType, PermutationIndex >::type IntColVectorType
 
typedef internal::plain_row_type< MatrixType, PermutationIndex >::type IntRowVectorType
 
typedef MatrixType_ MatrixType
 
using PermutationIndex = PermutationIndex_
 
typedef PermutationMatrix< RowsAtCompileTime, MaxRowsAtCompileTime, PermutationIndexPermutationPType
 
typedef PermutationMatrix< ColsAtCompileTime, MaxColsAtCompileTime, PermutationIndexPermutationQType
 
typedef MatrixType::PlainObject PlainObject
 
- Public Types inherited from Eigen::SolverBase< FullPivLU< MatrixType_, PermutationIndex_ > >
enum  
 
typedef std::conditional_t< NumTraits< Scalar >::IsComplex, CwiseUnaryOp< internal::scalar_conjugate_op< Scalar >, const ConstTransposeReturnType >, const ConstTransposeReturnTypeAdjointReturnType
 
typedef EigenBase< FullPivLU< MatrixType_, PermutationIndex_ > > Base
 
typedef Scalar CoeffReturnType
 
typedef Transpose< const FullPivLU< MatrixType_, PermutationIndex_ > > ConstTransposeReturnType
 
typedef internal::traits< FullPivLU< 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 >
FullPivLUcompute (const EigenBase< InputType > &matrix)
 
internal::traits< MatrixType >::Scalar determinant () const
 
Index dimensionOfKernel () const
 
 FullPivLU ()
 Default Constructor. 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...
 
 FullPivLU (Index rows, Index cols)
 Default Constructor with memory preallocation. 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
 
const PermutationPTypepermutationP () const
 
const PermutationQTypepermutationQ () const
 
Index rank () const
 
RealScalar rcond () const
 
MatrixType reconstructedMatrix () const
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
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_, PermutationIndex_ > >
const AdjointReturnType adjoint () const
 
FullPivLU< MatrixType_, PermutationIndex_ > & derived ()
 
const FullPivLU< MatrixType_, PermutationIndex_ > & derived () const
 
const Solve< FullPivLU< 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 computeInPlace ()
 
- Protected Member Functions inherited from Eigen::SolverBase< FullPivLU< MatrixType_, PermutationIndex_ > >
void _check_solve_assertion (const Rhs &b) const
 

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_, typename PermutationIndex_>
class Eigen::FullPivLU< MatrixType_, PermutationIndex_ >

LU decomposition of a matrix with complete 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 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 example, 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;
Matrix3f m
LU decomposition of a matrix with complete pivoting, and related features.
Definition: FullPivLU.h:64
cout<< "Here is the matrix m:"<< endl<< m<< endl;Eigen::FullPivLU< Matrix5x3 > lu(m)
Matrix< double, 5, 3 > Matrix5x3
Matrix< double, 5, 5 > Matrix5x5
@ Upper
Definition: Constants.h:213

Output:

Here is the matrix m:
   0.68  -0.605 -0.0452
 -0.211   -0.33   0.258
  0.566   0.536   -0.27
  0.597  -0.444  0.0268
  0.823   0.108   0.904
Here is, up to permutations, its LU decomposition matrix:
 0.904  0.823  0.108
-0.299  0.812  0.569
 -0.05  0.888   -1.1
0.0296  0.705  0.768
 0.285 -0.549 0.0436
Here is the L part:
     1      0      0      0      0
-0.299      1      0      0      0
 -0.05  0.888      1      0      0
0.0296  0.705  0.768      1      0
 0.285 -0.549 0.0436      0      1
Here is the U part:
0.904 0.823 0.108
    0 0.812 0.569
    0     0  -1.1
    0     0     0
    0     0     0
Let us now reconstruct the original matrix m:
   0.68  -0.605 -0.0452
 -0.211   -0.33   0.258
  0.566   0.536   -0.27
  0.597  -0.444  0.0268
  0.823   0.108   0.904

This class supports the inplace decomposition mechanism.

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

Definition at line 62 of file FullPivLU.h.

Member Typedef Documentation

◆ Base

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

Definition at line 67 of file FullPivLU.h.

◆ IntColVectorType

template<typename MatrixType_ , typename PermutationIndex_ >
typedef internal::plain_col_type<MatrixType, PermutationIndex>::type Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::IntColVectorType

Definition at line 77 of file FullPivLU.h.

◆ IntRowVectorType

template<typename MatrixType_ , typename PermutationIndex_ >
typedef internal::plain_row_type<MatrixType, PermutationIndex>::type Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::IntRowVectorType

Definition at line 76 of file FullPivLU.h.

◆ MatrixType

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

Definition at line 66 of file FullPivLU.h.

◆ PermutationIndex

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

Definition at line 75 of file FullPivLU.h.

◆ PermutationPType

template<typename MatrixType_ , typename PermutationIndex_ >
typedef PermutationMatrix<RowsAtCompileTime, MaxRowsAtCompileTime, PermutationIndex> Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::PermutationPType

Definition at line 79 of file FullPivLU.h.

◆ PermutationQType

template<typename MatrixType_ , typename PermutationIndex_ >
typedef PermutationMatrix<ColsAtCompileTime, MaxColsAtCompileTime, PermutationIndex> Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::PermutationQType

Definition at line 78 of file FullPivLU.h.

◆ PlainObject

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

Definition at line 80 of file FullPivLU.h.

Member Enumeration Documentation

◆ anonymous enum

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

Definition at line 71 of file FullPivLU.h.

71  {
72  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
73  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
74  };

Constructor & Destructor Documentation

◆ FullPivLU() [1/4]

template<typename MatrixType , typename PermutationIndex >
Eigen::FullPivLU< MatrixType, PermutationIndex >::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 442 of file FullPivLU.h.

444 {
445 }
bool m_usePrescribedThreshold
Definition: FullPivLU.h:438
bool m_isInitialized
Definition: FullPivLU.h:438

◆ FullPivLU() [2/4]

template<typename MatrixType , typename PermutationIndex >
Eigen::FullPivLU< MatrixType, PermutationIndex >::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 448 of file FullPivLU.h.

449  : m_lu(rows, cols),
450  m_p(rows),
451  m_q(cols),
454  m_isInitialized(false),
456 {
457 }
IntColVectorType m_rowsTranspositions
Definition: FullPivLU.h:432
MatrixType m_lu
Definition: FullPivLU.h:429
IntRowVectorType m_colsTranspositions
Definition: FullPivLU.h:433
EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: FullPivLU.h:411
PermutationPType m_p
Definition: FullPivLU.h:430
EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: FullPivLU.h:413
PermutationQType m_q
Definition: FullPivLU.h:431

◆ FullPivLU() [3/4]

template<typename MatrixType , typename PermutationIndex >
template<typename InputType >
Eigen::FullPivLU< MatrixType, PermutationIndex >::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 461 of file FullPivLU.h.

462  : m_lu(matrix.rows(), matrix.cols()),
463  m_p(matrix.rows()),
464  m_q(matrix.cols()),
465  m_rowsTranspositions(matrix.rows()),
466  m_colsTranspositions(matrix.cols()),
467  m_isInitialized(false),
469 {
470  compute(matrix.derived());
471 }
FullPivLU & compute(const EigenBase< InputType > &matrix)
Definition: FullPivLU.h:123

◆ FullPivLU() [4/4]

template<typename MatrixType , typename PermutationIndex >
template<typename InputType >
Eigen::FullPivLU< MatrixType, PermutationIndex >::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 475 of file FullPivLU.h.

476  : m_lu(matrix.derived()),
477  m_p(matrix.rows()),
478  m_q(matrix.cols()),
479  m_rowsTranspositions(matrix.rows()),
480  m_colsTranspositions(matrix.cols()),
481  m_isInitialized(false),
483 {
484  computeInPlace();
485 }
void computeInPlace()
Definition: FullPivLU.h:488

Member Function Documentation

◆ cols()

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

Definition at line 413 of file FullPivLU.h.

413 { return m_lu.cols(); }

◆ compute()

template<typename MatrixType_ , typename PermutationIndex_ >
template<typename InputType >
FullPivLU& Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::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 123 of file FullPivLU.h.

123  {
124  m_lu = matrix.derived();
125  computeInPlace();
126  return *this;
127  }

◆ computeInPlace()

template<typename MatrixType , typename PermutationIndex >
void Eigen::FullPivLU< MatrixType, PermutationIndex >::computeInPlace
protected

Definition at line 488 of file FullPivLU.h.

489 {
490  eigen_assert(m_lu.rows()<=NumTraits<PermutationIndex>::highest() && m_lu.cols()<=NumTraits<PermutationIndex>::highest());
491 
492  m_l1_norm = m_lu.cwiseAbs().colwise().sum().maxCoeff();
493 
494  const Index size = m_lu.diagonalSize();
495  const Index rows = m_lu.rows();
496  const Index cols = m_lu.cols();
497 
498  // will store the transpositions, before we accumulate them at the end.
499  // can't accumulate on-the-fly because that will be done in reverse order for the rows.
500  m_rowsTranspositions.resize(m_lu.rows());
501  m_colsTranspositions.resize(m_lu.cols());
502  Index number_of_transpositions = 0; // number of NONTRIVIAL transpositions, i.e. m_rowsTranspositions[i]!=i
503 
504  m_nonzero_pivots = size; // the generic case is that in which all pivots are nonzero (invertible case)
505  m_maxpivot = RealScalar(0);
506 
507  for(Index k = 0; k < size; ++k)
508  {
509  // First, we need to find the pivot.
510 
511  // biggest coefficient in the remaining bottom-right corner (starting at row k, col k)
512  Index row_of_biggest_in_corner, col_of_biggest_in_corner;
513  typedef internal::scalar_score_coeff_op<Scalar> Scoring;
514  typedef typename Scoring::result_type Score;
515  Score biggest_in_corner;
516  biggest_in_corner = m_lu.bottomRightCorner(rows-k, cols-k)
517  .unaryExpr(Scoring())
518  .maxCoeff(&row_of_biggest_in_corner, &col_of_biggest_in_corner);
519  row_of_biggest_in_corner += k; // correct the values! since they were computed in the corner,
520  col_of_biggest_in_corner += k; // need to add k to them.
521 
522  if(numext::is_exactly_zero(biggest_in_corner))
523  {
524  // before exiting, make sure to initialize the still uninitialized transpositions
525  // in a sane state without destroying what we already have.
526  m_nonzero_pivots = k;
527  for(Index i = k; i < size; ++i)
528  {
529  m_rowsTranspositions.coeffRef(i) = internal::convert_index<StorageIndex>(i);
530  m_colsTranspositions.coeffRef(i) = internal::convert_index<StorageIndex>(i);
531  }
532  break;
533  }
534 
535  RealScalar abs_pivot = internal::abs_knowing_score<Scalar>()(m_lu(row_of_biggest_in_corner, col_of_biggest_in_corner), biggest_in_corner);
536  if(abs_pivot > m_maxpivot) m_maxpivot = abs_pivot;
537 
538  // Now that we've found the pivot, we need to apply the row/col swaps to
539  // bring it to the location (k,k).
540 
541  m_rowsTranspositions.coeffRef(k) = internal::convert_index<StorageIndex>(row_of_biggest_in_corner);
542  m_colsTranspositions.coeffRef(k) = internal::convert_index<StorageIndex>(col_of_biggest_in_corner);
543  if(k != row_of_biggest_in_corner) {
544  m_lu.row(k).swap(m_lu.row(row_of_biggest_in_corner));
545  ++number_of_transpositions;
546  }
547  if(k != col_of_biggest_in_corner) {
548  m_lu.col(k).swap(m_lu.col(col_of_biggest_in_corner));
549  ++number_of_transpositions;
550  }
551 
552  // Now that the pivot is at the right location, we update the remaining
553  // bottom-right corner by Gaussian elimination.
554 
555  if(k<rows-1)
556  m_lu.col(k).tail(rows-k-1) /= m_lu.coeff(k,k);
557  if(k<size-1)
558  m_lu.block(k+1,k+1,rows-k-1,cols-k-1).noalias() -= m_lu.col(k).tail(rows-k-1) * m_lu.row(k).tail(cols-k-1);
559  }
560 
561  // the main loop is over, we still have to accumulate the transpositions to find the
562  // permutations P and Q
563 
565  for(Index k = size-1; k >= 0; --k)
567 
569  for(Index k = 0; k < size; ++k)
571 
572  m_det_pq = (number_of_transpositions%2) ? -1 : 1;
573 
574  m_isInitialized = true;
575 }
#define eigen_assert(x)
Definition: Macros.h:902
RealScalar m_l1_norm
Definition: FullPivLU.h:435
signed char m_det_pq
Definition: FullPivLU.h:437
Index m_nonzero_pivots
Definition: FullPivLU.h:434
RealScalar m_maxpivot
Definition: FullPivLU.h:436
Derived & applyTranspositionOnTheRight(Index i, Index j)
bool is_exactly_zero(const X &x)
Definition: Meta.h:475
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:41
EIGEN_CONSTEXPR Index size() const EIGEN_NOEXCEPT
Definition: EigenBase.h:69

◆ determinant()

template<typename MatrixType , typename PermutationIndex >
internal::traits< MatrixType >::Scalar Eigen::FullPivLU< 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
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 578 of file FullPivLU.h.

579 {
580  eigen_assert(m_isInitialized && "LU is not initialized.");
581  eigen_assert(m_lu.rows() == m_lu.cols() && "You can't take the determinant of a non-square matrix!");
582  return Scalar(m_det_pq) * Scalar(m_lu.diagonal().prod());
583 }
internal::traits< FullPivLU< MatrixType_, PermutationIndex_ > >::Scalar Scalar
Definition: SolverBase.h:75

◆ dimensionOfKernel()

template<typename MatrixType_ , typename PermutationIndex_ >
Index Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::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 350 of file FullPivLU.h.

351  {
352  eigen_assert(m_isInitialized && "LU is not initialized.");
353  return cols() - rank();
354  }
Index rank() const
Definition: FullPivLU.h:333

◆ image()

template<typename MatrixType_ , typename PermutationIndex_ >
const internal::image_retval<FullPivLU> Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::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:

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;
Matrix< double, 3, 3 > Matrix3d
3×3 matrix of type double.
Definition: Matrix.h:502

Output:

Here is the matrix m:
1 1 0
1 3 2
0 1 1
Notice that the middle column is the sum of the two others, so the columns are linearly dependent.
Here is a matrix whose columns have the same span but are linearly independent:
1 1
3 1
1 0
See also
kernel()

Definition at line 219 of file FullPivLU.h.

220  {
221  eigen_assert(m_isInitialized && "LU is not initialized.");
222  return internal::image_retval<FullPivLU>(*this, originalMatrix);
223  }

◆ inverse()

template<typename MatrixType_ , typename PermutationIndex_ >
const Inverse<FullPivLU> Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::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 401 of file FullPivLU.h.

402  {
403  eigen_assert(m_isInitialized && "LU is not initialized.");
404  eigen_assert(m_lu.rows() == m_lu.cols() && "You can't take the inverse of a non-square matrix!");
405  return Inverse<FullPivLU>(*this);
406  }

◆ isInjective()

template<typename MatrixType_ , typename PermutationIndex_ >
bool Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::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 363 of file FullPivLU.h.

364  {
365  eigen_assert(m_isInitialized && "LU is not initialized.");
366  return rank() == cols();
367  }

◆ isInvertible()

template<typename MatrixType_ , typename PermutationIndex_ >
bool Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::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 388 of file FullPivLU.h.

389  {
390  eigen_assert(m_isInitialized && "LU is not initialized.");
391  return isInjective() && (m_lu.rows() == m_lu.cols());
392  }
bool isInjective() const
Definition: FullPivLU.h:363

◆ isSurjective()

template<typename MatrixType_ , typename PermutationIndex_ >
bool Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::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 376 of file FullPivLU.h.

377  {
378  eigen_assert(m_isInitialized && "LU is not initialized.");
379  return rank() == rows();
380  }

◆ kernel()

template<typename MatrixType_ , typename PermutationIndex_ >
const internal::kernel_retval<FullPivLU> Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::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:

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;
static const RandomReturnType Random()
Definition: Random.h:114
Matrix< float, Dynamic, Dynamic > MatrixXf
Dynamic×Dynamic matrix of type float.
Definition: Matrix.h:501

Output:

Here is the matrix m:
   0.68   0.597   -0.33   0.108   -0.27
 -0.211   0.823   0.536 -0.0452  0.0268
  0.566  -0.605  -0.444   0.258   0.904
Here is a matrix whose columns form a basis of the kernel of m:
 -0.219   0.763
0.00335  -0.447
      0       1
      1       0
 -0.145  -0.285
By definition of the kernel, m*ker is zero:
 7.45e-09  1.49e-08
-1.86e-09 -4.05e-08
        0 -2.98e-08
See also
image()

Definition at line 193 of file FullPivLU.h.

194  {
195  eigen_assert(m_isInitialized && "LU is not initialized.");
196  return internal::kernel_retval<FullPivLU>(*this);
197  }

◆ matrixLU()

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

136  {
137  eigen_assert(m_isInitialized && "LU is not initialized.");
138  return m_lu;
139  }

◆ maxPivot()

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

Definition at line 157 of file FullPivLU.h.

157 { return m_maxpivot; }

◆ nonzeroPivots()

template<typename MatrixType_ , typename PermutationIndex_ >
Index Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::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 148 of file FullPivLU.h.

149  {
150  eigen_assert(m_isInitialized && "LU is not initialized.");
151  return m_nonzero_pivots;
152  }

◆ permutationP()

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

Definition at line 163 of file FullPivLU.h.

164  {
165  eigen_assert(m_isInitialized && "LU is not initialized.");
166  return m_p;
167  }

◆ permutationQ()

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

Definition at line 173 of file FullPivLU.h.

174  {
175  eigen_assert(m_isInitialized && "LU is not initialized.");
176  return m_q;
177  }

◆ rank()

template<typename MatrixType_ , typename PermutationIndex_ >
Index Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::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 333 of file FullPivLU.h.

334  {
335  using std::abs;
336  eigen_assert(m_isInitialized && "LU is not initialized.");
337  RealScalar premultiplied_threshold = abs(m_maxpivot) * threshold();
338  Index result = 0;
339  for(Index i = 0; i < m_nonzero_pivots; ++i)
340  result += (abs(m_lu.coeff(i,i)) > premultiplied_threshold);
341  return result;
342  }
const AbsReturnType abs() const
RealScalar threshold() const
Definition: FullPivLU.h:318
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)

◆ rcond()

template<typename MatrixType_ , typename PermutationIndex_ >
RealScalar Eigen::FullPivLU< 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 253 of file FullPivLU.h.

254  {
255  eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
257  }
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::FullPivLU< MatrixType, PermutationIndex >::reconstructedMatrix
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 589 of file FullPivLU.h.

590 {
591  eigen_assert(m_isInitialized && "LU is not initialized.");
592  const Index smalldim = (std::min)(m_lu.rows(), m_lu.cols());
593  // LU
594  MatrixType res(m_lu.rows(),m_lu.cols());
595  // FIXME the .toDenseMatrix() should not be needed...
596  res = m_lu.leftCols(smalldim)
597  .template triangularView<UnitLower>().toDenseMatrix()
598  * m_lu.topRows(smalldim)
599  .template triangularView<Upper>().toDenseMatrix();
600 
601  // P^{-1}(LU)
602  res = m_p.inverse() * res;
603 
604  // (P^{-1}LU)Q^{-1}
605  res = res * m_q.inverse();
606 
607  return res;
608 }
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
Matrix< float, 1, Dynamic > MatrixType
InverseReturnType inverse() const
bfloat16() min(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:684

◆ rows()

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

Definition at line 411 of file FullPivLU.h.

411 { return m_lu.rows(); }

◆ setThreshold() [1/2]

template<typename MatrixType_ , typename PermutationIndex_ >
FullPivLU& Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::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 293 of file FullPivLU.h.

294  {
297  return *this;
298  }
RealScalar m_prescribedThreshold
Definition: FullPivLU.h:436

◆ setThreshold() [2/2]

template<typename MatrixType_ , typename PermutationIndex_ >
FullPivLU& Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::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);
@ Default
Definition: Constants.h:364

See the documentation of setThreshold(const RealScalar&).

Definition at line 308 of file FullPivLU.h.

309  {
310  m_usePrescribedThreshold = false;
311  return *this;
312  }

◆ solve()

template<typename MatrixType_ , typename PermutationIndex_ >
template<typename Rhs >
const Solve<FullPivLU, Rhs> Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::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.

This method just tries to find as good a solution as possible. If you want to check whether a solution exists or if it is accurate, just call this function to get a result and then compute the error of this result, or use MatrixBase::isApprox() directly, for instance like this:

bool a_solution_exists = (A*result).isApprox(b, precision);
Array< int, 3, 1 > b
MatrixXcf A
bool isApprox(const Scalar &x, const Scalar &y, const typename NumTraits< Scalar >::Real &precision=NumTraits< Scalar >::dummy_precision())

This method avoids dividing by zero, so that the non-existence of a solution doesn't by itself mean that you'll get inf or nan values.

If there exists more than one solution, this method will arbitrarily choose one. If you need a complete analysis of the space of solutions, take the one solution obtained by this method and add to it elements of the kernel, as determined by kernel().

Example:

Matrix<float,2,3> m = Matrix<float,2,3>::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;
Matrix< float, 2, 2 > Matrix2f
2×2 matrix of type float.
Definition: Matrix.h:501
const Scalar & y

Output:

Here is the matrix m:
  0.68  0.566  0.823
-0.211  0.597 -0.605
Here is the matrix y:
 -0.33 -0.444
 0.536  0.108
Here is a solution x to the equation mx=y:
     0      0
 0.291 -0.216
  -0.6 -0.391
See also
TriangularView::solve(), kernel(), inverse()

◆ threshold()

template<typename MatrixType_ , typename PermutationIndex_ >
RealScalar Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::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 318 of file FullPivLU.h.

319  {
322  // this formula comes from experimenting (see "LU precision tuning" thread on the list)
323  // and turns out to be identical to Higham's formula used already in LDLt.
324  : NumTraits<Scalar>::epsilon() * RealScalar(m_lu.diagonalSize());
325  }

Member Data Documentation

◆ m_colsTranspositions

template<typename MatrixType_ , typename PermutationIndex_ >
IntRowVectorType Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_colsTranspositions
protected

Definition at line 433 of file FullPivLU.h.

◆ m_det_pq

template<typename MatrixType_ , typename PermutationIndex_ >
signed char Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_det_pq
protected

Definition at line 437 of file FullPivLU.h.

◆ m_isInitialized

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

Definition at line 438 of file FullPivLU.h.

◆ m_l1_norm

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

Definition at line 435 of file FullPivLU.h.

◆ m_lu

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

Definition at line 429 of file FullPivLU.h.

◆ m_maxpivot

template<typename MatrixType_ , typename PermutationIndex_ >
RealScalar Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_maxpivot
protected

Definition at line 436 of file FullPivLU.h.

◆ m_nonzero_pivots

template<typename MatrixType_ , typename PermutationIndex_ >
Index Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_nonzero_pivots
protected

Definition at line 434 of file FullPivLU.h.

◆ m_p

template<typename MatrixType_ , typename PermutationIndex_ >
PermutationPType Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_p
protected

Definition at line 430 of file FullPivLU.h.

◆ m_prescribedThreshold

template<typename MatrixType_ , typename PermutationIndex_ >
RealScalar Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_prescribedThreshold
protected

Definition at line 436 of file FullPivLU.h.

◆ m_q

template<typename MatrixType_ , typename PermutationIndex_ >
PermutationQType Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_q
protected

Definition at line 431 of file FullPivLU.h.

◆ m_rowsTranspositions

template<typename MatrixType_ , typename PermutationIndex_ >
IntColVectorType Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_rowsTranspositions
protected

Definition at line 432 of file FullPivLU.h.

◆ m_usePrescribedThreshold

template<typename MatrixType_ , typename PermutationIndex_ >
bool Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_usePrescribedThreshold
protected

Definition at line 438 of file FullPivLU.h.


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