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

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

#include <FullPivLU.h>

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

Public Types

enum  { MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime , MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime }
 
typedef MatrixType_ MatrixType
 
typedef SolverBase< FullPivLUBase
 
using PermutationIndex = PermutationIndex_
 
typedef internal::plain_row_type< MatrixType, PermutationIndex >::type IntRowVectorType
 
typedef internal::plain_col_type< MatrixType, PermutationIndex >::type IntColVectorType
 
typedef PermutationMatrix< ColsAtCompileTime, MaxColsAtCompileTime, PermutationIndexPermutationQType
 
typedef PermutationMatrix< RowsAtCompileTime, MaxRowsAtCompileTime, PermutationIndexPermutationPType
 
typedef MatrixType::PlainObject PlainObject
 
- Public Types inherited from Eigen::SolverBase< FullPivLU< MatrixType_, PermutationIndex_ > >
enum  
 
typedef EigenBase< FullPivLU< MatrixType_, PermutationIndex_ > > Base
 
typedef internal::traits< FullPivLU< MatrixType_, PermutationIndex_ > >::Scalar Scalar
 
typedef Scalar CoeffReturnType
 
typedef Transpose< const FullPivLU< MatrixType_, PermutationIndex_ > > ConstTransposeReturnType
 
typedef std::conditional_t< NumTraits< Scalar >::IsComplex, CwiseUnaryOp< internal::scalar_conjugate_op< Scalar >, const ConstTransposeReturnType >, const ConstTransposeReturnTypeAdjointReturnType
 
- 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

 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...
 
template<typename InputType >
FullPivLUcompute (const EigenBase< InputType > &matrix)
 
const MatrixTypematrixLU () const
 
Index nonzeroPivots () const
 
RealScalar maxPivot () const
 
EIGEN_DEVICE_FUNC const PermutationPTypepermutationP () const
 
const PermutationQTypepermutationQ () const
 
const internal::kernel_retval< FullPivLUkernel () const
 
const internal::image_retval< FullPivLUimage (const MatrixType &originalMatrix) const
 
RealScalar rcond () const
 
internal::traits< MatrixType >::Scalar determinant () const
 
FullPivLUsetThreshold (const RealScalar &threshold)
 
FullPivLUsetThreshold (Default_t)
 
RealScalar threshold () const
 
Index rank () const
 
Index dimensionOfKernel () const
 
bool isInjective () const
 
bool isSurjective () const
 
bool isInvertible () const
 
const Inverse< FullPivLUinverse () const
 
MatrixType reconstructedMatrix () const
 
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
template<typename RhsType , typename DstType >
void _solve_impl (const RhsType &rhs, DstType &dst) const
 
template<bool Conjugate, typename RhsType , typename DstType >
void _solve_impl_transposed (const RhsType &rhs, DstType &dst) const
 
- Public Member Functions inherited from Eigen::SolverBase< FullPivLU< MatrixType_, PermutationIndex_ > >
 SolverBase ()
 
 ~SolverBase ()
 
const Solve< FullPivLU< MatrixType_, PermutationIndex_ >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
const ConstTransposeReturnType transpose () const
 
const AdjointReturnType adjoint () const
 
constexpr EIGEN_DEVICE_FUNC FullPivLU< MatrixType_, PermutationIndex_ > & derived ()
 
constexpr EIGEN_DEVICE_FUNC const FullPivLU< MatrixType_, PermutationIndex_ > & derived () const
 
- Public Member Functions inherited from Eigen::EigenBase< Derived >
constexpr EIGEN_DEVICE_FUNC Derived & derived ()
 
constexpr EIGEN_DEVICE_FUNC const Derived & derived () const
 
EIGEN_DEVICE_FUNC Derived & const_cast_derived () const
 
EIGEN_DEVICE_FUNC const Derived & const_derived () const
 
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index size () const EIGEN_NOEXCEPT
 
template<typename Dest >
EIGEN_DEVICE_FUNC void evalTo (Dest &dst) const
 
template<typename Dest >
EIGEN_DEVICE_FUNC void addTo (Dest &dst) const
 
template<typename Dest >
EIGEN_DEVICE_FUNC void subTo (Dest &dst) const
 
template<typename Dest >
EIGEN_DEVICE_FUNC void applyThisOnTheRight (Dest &dst) const
 
template<typename Dest >
EIGEN_DEVICE_FUNC void applyThisOnTheLeft (Dest &dst) const
 
template<typename Device >
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE DeviceWrapper< Derived, Device > device (Device &device)
 
template<typename Device >
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE DeviceWrapper< const Derived, Device > device (Device &device) 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

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

Friends

class SolverBase< FullPivLU >
 

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;
LU decomposition of a matrix with complete pivoting, and related features.
Definition: FullPivLU.h:63
cout<< "Here is the matrix m:"<< endl<< m<< endl;Eigen::FullPivLU< Matrix5x3 > lu(m)
Matrix< double, 5, 3 > Matrix5x3
Definition: class_FullPivLU.cpp:1
Matrix< double, 5, 5 > Matrix5x5
Definition: class_FullPivLU.cpp:2
@ Upper
Definition: Constants.h:213
int * m
Definition: level2_cplx_impl.h:294

Output:

This class supports the inplace decomposition mechanism.

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

Member Typedef Documentation

◆ Base

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

◆ IntColVectorType

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

◆ IntRowVectorType

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

◆ MatrixType

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

◆ PermutationIndex

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

◆ PermutationPType

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

◆ PermutationQType

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

◆ PlainObject

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

Member Enumeration Documentation

◆ anonymous enum

template<typename MatrixType_ , typename PermutationIndex_ >
anonymous enum
Enumerator
MaxRowsAtCompileTime 
MaxColsAtCompileTime 
70  {
71  MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
72  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
73  };
@ MaxRowsAtCompileTime
Definition: FullPivLU.h:71
@ MaxColsAtCompileTime
Definition: FullPivLU.h:72

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

420 : m_isInitialized(false), m_usePrescribedThreshold(false) {}
bool m_usePrescribedThreshold
Definition: FullPivLU.h:416
bool m_isInitialized
Definition: FullPivLU.h:416

◆ 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()
424  : m_lu(rows, cols),
425  m_p(rows),
426  m_q(cols),
429  m_isInitialized(false),
430  m_usePrescribedThreshold(false) {}
IntColVectorType m_rowsTranspositions
Definition: FullPivLU.h:410
MatrixType m_lu
Definition: FullPivLU.h:407
IntRowVectorType m_colsTranspositions
Definition: FullPivLU.h:411
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: FullPivLU.h:391
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: FullPivLU.h:392
PermutationPType m_p
Definition: FullPivLU.h:408
PermutationQType m_q
Definition: FullPivLU.h:409

◆ 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.
435  : m_lu(matrix.rows(), matrix.cols()),
436  m_p(matrix.rows()),
437  m_q(matrix.cols()),
440  m_isInitialized(false),
441  m_usePrescribedThreshold(false) {
442  compute(matrix.derived());
443 }
FullPivLU & compute(const EigenBase< InputType > &matrix)
Definition: FullPivLU.h:123
Eigen::Map< Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor >, 0, Eigen::OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Definition: common.h:85

References Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::compute(), and matrix().

◆ 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&)
448  : m_lu(matrix.derived()),
449  m_p(matrix.rows()),
450  m_q(matrix.cols()),
453  m_isInitialized(false),
454  m_usePrescribedThreshold(false) {
455  computeInPlace();
456 }
void computeInPlace()
Definition: FullPivLU.h:459

References Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::computeInPlace().

Member Function Documentation

◆ _solve_impl()

template<typename MatrixType_ , typename PermutationIndex_ >
template<typename RhsType , typename DstType >
void Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::_solve_impl ( const RhsType &  rhs,
DstType &  dst 
) const
693  {
694  /* The decomposition PAQ = LU can be rewritten as A = P^{-1} L U Q^{-1}.
695  * So we proceed as follows:
696  * Step 1: compute c = P * rhs.
697  * Step 2: replace c by the solution x to Lx = c. Exists because L is invertible.
698  * Step 3: replace c by the solution x to Ux = c. May or may not exist.
699  * Step 4: result = Q * c;
700  */
701 
702  const Index rows = this->rows(), cols = this->cols(), nonzero_pivots = this->rank();
703  const Index smalldim = (std::min)(rows, cols);
704 
705  if (nonzero_pivots == 0) {
706  dst.setZero();
707  return;
708  }
709 
710  typename RhsType::PlainObject c(rhs.rows(), rhs.cols());
711 
712  // Step 1
713  c = permutationP() * rhs;
714 
715  // Step 2
716  m_lu.topLeftCorner(smalldim, smalldim).template triangularView<UnitLower>().solveInPlace(c.topRows(smalldim));
717  if (rows > cols) c.bottomRows(rows - cols) -= m_lu.bottomRows(rows - cols) * c.topRows(cols);
718 
719  // Step 3
720  m_lu.topLeftCorner(nonzero_pivots, nonzero_pivots)
721  .template triangularView<Upper>()
722  .solveInPlace(c.topRows(nonzero_pivots));
723 
724  // Step 4
725  for (Index i = 0; i < nonzero_pivots; ++i) dst.row(permutationQ().indices().coeff(i)) = c.row(i);
726  for (Index i = nonzero_pivots; i < m_lu.cols(); ++i) dst.row(permutationQ().indices().coeff(i)).setZero();
727 }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Index rank() const
Definition: FullPivLU.h:321
EIGEN_DEVICE_FUNC const PermutationPType & permutationP() const
Definition: FullPivLU.h:161
const PermutationQType & permutationQ() const
Definition: FullPivLU.h:170
#define min(a, b)
Definition: datatypes.h:22
int c
Definition: calibrate.py:100
Eigen::Index Index
The interface type of indices.
Definition: EigenBase.h:43

References calibrate::c, cols, i, min, and rows.

◆ _solve_impl_transposed()

template<typename MatrixType_ , typename PermutationIndex_ >
template<bool Conjugate, typename RhsType , typename DstType >
void Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::_solve_impl_transposed ( const RhsType &  rhs,
DstType &  dst 
) const
731  {
732  /* The decomposition PAQ = LU can be rewritten as A = P^{-1} L U Q^{-1},
733  * and since permutations are real and unitary, we can write this
734  * as A^T = Q U^T L^T P,
735  * So we proceed as follows:
736  * Step 1: compute c = Q^T rhs.
737  * Step 2: replace c by the solution x to U^T x = c. May or may not exist.
738  * Step 3: replace c by the solution x to L^T x = c.
739  * Step 4: result = P^T c.
740  * If Conjugate is true, replace "^T" by "^*" above.
741  */
742 
743  const Index rows = this->rows(), cols = this->cols(), nonzero_pivots = this->rank();
744  const Index smalldim = (std::min)(rows, cols);
745 
746  if (nonzero_pivots == 0) {
747  dst.setZero();
748  return;
749  }
750 
751  typename RhsType::PlainObject c(rhs.rows(), rhs.cols());
752 
753  // Step 1
754  c = permutationQ().inverse() * rhs;
755 
756  // Step 2
757  m_lu.topLeftCorner(nonzero_pivots, nonzero_pivots)
758  .template triangularView<Upper>()
759  .transpose()
760  .template conjugateIf<Conjugate>()
761  .solveInPlace(c.topRows(nonzero_pivots));
762 
763  // Step 3
764  m_lu.topLeftCorner(smalldim, smalldim)
765  .template triangularView<UnitLower>()
766  .transpose()
767  .template conjugateIf<Conjugate>()
768  .solveInPlace(c.topRows(smalldim));
769 
770  // Step 4
771  PermutationPType invp = permutationP().inverse().eval();
772  for (Index i = 0; i < smalldim; ++i) dst.row(invp.indices().coeff(i)) = c.row(i);
773  for (Index i = smalldim; i < rows; ++i) dst.row(invp.indices().coeff(i)).setZero();
774 }
PermutationMatrix< RowsAtCompileTime, MaxRowsAtCompileTime, PermutationIndex > PermutationPType
Definition: FullPivLU.h:78
InverseReturnType inverse() const
Definition: PermutationMatrix.h:172

References calibrate::c, cols, i, Eigen::PermutationMatrix< SizeAtCompileTime, MaxSizeAtCompileTime, StorageIndex_ >::indices(), Eigen::PermutationBase< Derived >::inverse(), min, and rows.

◆ 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
123  {
124  m_lu = matrix.derived();
125  computeInPlace();
126  return *this;
127  }

References Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::computeInPlace(), Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_lu, and matrix().

Referenced by ctms_decompositions(), Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::FullPivLU(), Eigen::internal::idrs(), and Eigen::internal::idrstabl().

◆ computeInPlace()

template<typename MatrixType , typename PermutationIndex >
void Eigen::FullPivLU< MatrixType, PermutationIndex >::computeInPlace
protected
459  {
460  eigen_assert(m_lu.rows() <= NumTraits<PermutationIndex>::highest() &&
461  m_lu.cols() <= NumTraits<PermutationIndex>::highest());
462 
463  m_l1_norm = m_lu.cwiseAbs().colwise().sum().maxCoeff();
464 
465  const Index size = m_lu.diagonalSize();
466  const Index rows = m_lu.rows();
467  const Index cols = m_lu.cols();
468 
469  // will store the transpositions, before we accumulate them at the end.
470  // can't accumulate on-the-fly because that will be done in reverse order for the rows.
471  m_rowsTranspositions.resize(m_lu.rows());
472  m_colsTranspositions.resize(m_lu.cols());
473  Index number_of_transpositions = 0; // number of NONTRIVIAL transpositions, i.e. m_rowsTranspositions[i]!=i
474 
475  m_nonzero_pivots = size; // the generic case is that in which all pivots are nonzero (invertible case)
476  m_maxpivot = RealScalar(0);
477 
478  for (Index k = 0; k < size; ++k) {
479  // First, we need to find the pivot.
480 
481  // biggest coefficient in the remaining bottom-right corner (starting at row k, col k)
482  Index row_of_biggest_in_corner, col_of_biggest_in_corner;
483  typedef internal::scalar_score_coeff_op<Scalar> Scoring;
484  typedef typename Scoring::result_type Score;
485  Score biggest_in_corner;
486  biggest_in_corner = m_lu.bottomRightCorner(rows - k, cols - k)
487  .unaryExpr(Scoring())
488  .maxCoeff(&row_of_biggest_in_corner, &col_of_biggest_in_corner);
489  row_of_biggest_in_corner += k; // correct the values! since they were computed in the corner,
490  col_of_biggest_in_corner += k; // need to add k to them.
491 
492  if (numext::is_exactly_zero(biggest_in_corner)) {
493  // before exiting, make sure to initialize the still uninitialized transpositions
494  // in a sane state without destroying what we already have.
496  for (Index i = k; i < size; ++i) {
497  m_rowsTranspositions.coeffRef(i) = internal::convert_index<StorageIndex>(i);
498  m_colsTranspositions.coeffRef(i) = internal::convert_index<StorageIndex>(i);
499  }
500  break;
501  }
502 
503  RealScalar abs_pivot = internal::abs_knowing_score<Scalar>()(
504  m_lu(row_of_biggest_in_corner, col_of_biggest_in_corner), biggest_in_corner);
505  if (abs_pivot > m_maxpivot) m_maxpivot = abs_pivot;
506 
507  // Now that we've found the pivot, we need to apply the row/col swaps to
508  // bring it to the location (k,k).
509 
510  m_rowsTranspositions.coeffRef(k) = internal::convert_index<StorageIndex>(row_of_biggest_in_corner);
511  m_colsTranspositions.coeffRef(k) = internal::convert_index<StorageIndex>(col_of_biggest_in_corner);
512  if (k != row_of_biggest_in_corner) {
513  m_lu.row(k).swap(m_lu.row(row_of_biggest_in_corner));
514  ++number_of_transpositions;
515  }
516  if (k != col_of_biggest_in_corner) {
517  m_lu.col(k).swap(m_lu.col(col_of_biggest_in_corner));
518  ++number_of_transpositions;
519  }
520 
521  // Now that the pivot is at the right location, we update the remaining
522  // bottom-right corner by Gaussian elimination.
523 
524  if (k < rows - 1) m_lu.col(k).tail(rows - k - 1) /= m_lu.coeff(k, k);
525  if (k < size - 1)
526  m_lu.block(k + 1, k + 1, rows - k - 1, cols - k - 1).noalias() -=
527  m_lu.col(k).tail(rows - k - 1) * m_lu.row(k).tail(cols - k - 1);
528  }
529 
530  // the main loop is over, we still have to accumulate the transpositions to find the
531  // permutations P and Q
532 
534  for (Index k = size - 1; k >= 0; --k) m_p.applyTranspositionOnTheRight(k, m_rowsTranspositions.coeff(k));
535 
538 
539  m_det_pq = (number_of_transpositions % 2) ? -1 : 1;
540 
541  m_isInitialized = true;
542 }
#define eigen_assert(x)
Definition: Macros.h:910
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:46
RealScalar m_l1_norm
Definition: FullPivLU.h:413
signed char m_det_pq
Definition: FullPivLU.h:415
Index m_nonzero_pivots
Definition: FullPivLU.h:412
RealScalar m_maxpivot
Definition: FullPivLU.h:414
Derived & applyTranspositionOnTheRight(Index i, Index j)
Definition: PermutationMatrix.h:162
void setIdentity()
Definition: PermutationMatrix.h:122
char char char int int * k
Definition: level2_impl.h:374
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bool is_exactly_zero(const X &x)
Definition: Meta.h:592
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index size() const EIGEN_NOEXCEPT
Definition: EigenBase.h:64

References cols, eigen_assert, i, Eigen::numext::is_exactly_zero(), k, rows, and size.

Referenced by Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::compute(), and Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::FullPivLU().

◆ 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()
545  {
546  eigen_assert(m_isInitialized && "LU is not initialized.");
547  eigen_assert(m_lu.rows() == m_lu.cols() && "You can't take the determinant of a non-square matrix!");
548  return Scalar(m_det_pq) * Scalar(m_lu.diagonal().prod());
549 }
internal::traits< FullPivLU< MatrixType_, PermutationIndex_ > >::Scalar Scalar
Definition: SolverBase.h:75

References eigen_assert.

◆ 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&).
336  {
337  eigen_assert(m_isInitialized && "LU is not initialized.");
338  return cols() - rank();
339  }

References Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::cols(), eigen_assert, Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_isInitialized, and Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::rank().

◆ 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:

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()
213  {
214  eigen_assert(m_isInitialized && "LU is not initialized.");
215  return internal::image_retval<FullPivLU>(*this, originalMatrix);
216  }

References eigen_assert, and Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_isInitialized.

Referenced by main().

◆ 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()
383  {
384  eigen_assert(m_isInitialized && "LU is not initialized.");
385  eigen_assert(m_lu.rows() == m_lu.cols() && "You can't take the inverse of a non-square matrix!");
386  return Inverse<FullPivLU>(*this);
387  }

References eigen_assert, Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_isInitialized, and Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_lu.

◆ 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&).
348  {
349  eigen_assert(m_isInitialized && "LU is not initialized.");
350  return rank() == cols();
351  }

References Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::cols(), eigen_assert, Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_isInitialized, and Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::rank().

Referenced by Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::isInvertible().

◆ 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&).
371  {
372  eigen_assert(m_isInitialized && "LU is not initialized.");
373  return isInjective() && (m_lu.rows() == m_lu.cols());
374  }
bool isInjective() const
Definition: FullPivLU.h:348

References eigen_assert, Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::isInjective(), Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_isInitialized, and Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_lu.

Referenced by inverse_general_4x4().

◆ 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&).
360  {
361  eigen_assert(m_isInitialized && "LU is not initialized.");
362  return rank() == rows();
363  }

References eigen_assert, Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_isInitialized, Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::rank(), and Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::rows().

◆ 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:

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()
189  {
190  eigen_assert(m_isInitialized && "LU is not initialized.");
191  return internal::kernel_retval<FullPivLU>(*this);
192  }

References eigen_assert, and Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_isInitialized.

Referenced by main().

◆ 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()
135  {
136  eigen_assert(m_isInitialized && "LU is not initialized.");
137  return m_lu;
138  }

References eigen_assert, Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_isInitialized, and Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_lu.

◆ 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.
155 { return m_maxpivot; }

References Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::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()
147  {
148  eigen_assert(m_isInitialized && "LU is not initialized.");
149  return m_nonzero_pivots;
150  }

References eigen_assert, Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_isInitialized, and Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_nonzero_pivots.

◆ permutationP()

template<typename MatrixType_ , typename PermutationIndex_ >
EIGEN_DEVICE_FUNC const PermutationPType& Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::permutationP ( ) const
inline
Returns
the permutation matrix P
See also
permutationQ()
161  {
162  eigen_assert(m_isInitialized && "LU is not initialized.");
163  return m_p;
164  }

References eigen_assert, Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_isInitialized, and Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_p.

◆ permutationQ()

template<typename MatrixType_ , typename PermutationIndex_ >
const PermutationQType& Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::permutationQ ( ) const
inline
Returns
the permutation matrix Q
See also
permutationP()
170  {
171  eigen_assert(m_isInitialized && "LU is not initialized.");
172  return m_q;
173  }

References eigen_assert, Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_isInitialized, and Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_q.

◆ 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&).
321  {
322  using std::abs;
323  eigen_assert(m_isInitialized && "LU is not initialized.");
324  RealScalar premultiplied_threshold = abs(m_maxpivot) * threshold();
325  Index result = 0;
326  for (Index i = 0; i < m_nonzero_pivots; ++i) result += (abs(m_lu.coeff(i, i)) > premultiplied_threshold);
327  return result;
328  }
AnnoyingScalar abs(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:135
RealScalar threshold() const
Definition: FullPivLU.h:307

References abs(), eigen_assert, i, Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_isInitialized, Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_lu, Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_maxpivot, Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_nonzero_pivots, and Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::threshold().

Referenced by Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::dimensionOfKernel(), Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::isInjective(), Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::isSurjective(), and main().

◆ 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.
245  {
246  eigen_assert(m_isInitialized && "PartialPivLU is not initialized.");
248  }
Decomposition::RealScalar rcond_estimate_helper(typename Decomposition::RealScalar matrix_norm, const Decomposition &dec)
Reciprocal condition number estimator.
Definition: ConditionEstimator.h:157

References eigen_assert, Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_isInitialized, Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_l1_norm, and Eigen::internal::rcond_estimate_helper().

◆ 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.
555  {
556  eigen_assert(m_isInitialized && "LU is not initialized.");
557  const Index smalldim = (std::min)(m_lu.rows(), m_lu.cols());
558  // LU
559  MatrixType res(m_lu.rows(), m_lu.cols());
560  // FIXME the .toDenseMatrix() should not be needed...
561  res = m_lu.leftCols(smalldim).template triangularView<UnitLower>().toDenseMatrix() *
562  m_lu.topRows(smalldim).template triangularView<Upper>().toDenseMatrix();
563 
564  // P^{-1}(LU)
565  res = m_p.inverse() * res;
566 
567  // (P^{-1}LU)Q^{-1}
568  res = res * m_q.inverse();
569 
570  return res;
571 }
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
Definition: PartialRedux_count.cpp:3
MatrixXf MatrixType
Definition: benchmark-blocking-sizes.cpp:52

References eigen_assert, min, and res.

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

284  {
287  return *this;
288  }
RealScalar m_prescribedThreshold
Definition: FullPivLU.h:414

References Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_prescribedThreshold, Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_usePrescribedThreshold, and Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::threshold().

◆ 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:361

See the documentation of setThreshold(const RealScalar&).

298  {
299  m_usePrescribedThreshold = false;
300  return *this;
301  }

References Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_usePrescribedThreshold.

◆ 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&).

307  {
310  // this formula comes from experimenting (see "LU precision tuning" thread on the
311  // list) and turns out to be identical to Higham's formula used already in LDLt.
312  : NumTraits<Scalar>::epsilon() * RealScalar(m_lu.diagonalSize());
313  }
double epsilon
Definition: osc_ring_sarah_asymptotics.h:43

References eigen_assert, oomph::SarahBL::epsilon, Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_isInitialized, Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_lu, Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_prescribedThreshold, and Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::m_usePrescribedThreshold.

Referenced by Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::rank(), and Eigen::FullPivLU< MatrixType_, PermutationIndex_ >::setThreshold().

Friends And Related Function Documentation

◆ SolverBase< FullPivLU >

template<typename MatrixType_ , typename PermutationIndex_ >
friend class SolverBase< FullPivLU >
friend

Member Data Documentation

◆ m_colsTranspositions

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

◆ m_det_pq

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

◆ m_isInitialized

◆ m_l1_norm

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

◆ m_lu

◆ m_maxpivot

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

◆ m_nonzero_pivots

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

◆ m_p

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

◆ m_prescribedThreshold

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

◆ m_q

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

◆ m_rowsTranspositions

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

◆ m_usePrescribedThreshold

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

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