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

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

#include <FullPivLU.h>

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

Public Types

enum  {
  RowsAtCompileTime = MatrixType::RowsAtCompileTime, ColsAtCompileTime = MatrixType::ColsAtCompileTime, Options = MatrixType::Options, MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
}
 
typedef _MatrixType MatrixType
 
typedef MatrixType::Scalar Scalar
 
typedef NumTraits< typename MatrixType::Scalar >::Real RealScalar
 
typedef internal::traits< MatrixType >::StorageKind StorageKind
 
typedef MatrixType::Index Index
 
typedef internal::plain_row_type< MatrixType, Index >::type IntRowVectorType
 
typedef internal::plain_col_type< MatrixType, Index >::type IntColVectorType
 
typedef PermutationMatrix< ColsAtCompileTime, MaxColsAtCompileTime > PermutationQType
 
typedef PermutationMatrix< RowsAtCompileTime, MaxRowsAtCompileTime > PermutationPType
 

Public Member Functions

 FullPivLU ()
 Default Constructor. More...
 
 FullPivLU (Index rows, Index cols)
 Default Constructor with memory preallocation. More...
 
 FullPivLU (const MatrixType &matrix)
 
FullPivLUcompute (const MatrixType &matrix)
 
const MatrixType & matrixLU () const
 
Index nonzeroPivots () const
 
RealScalar maxPivot () const
 
const PermutationPTypepermutationP () const
 
const PermutationQTypepermutationQ () const
 
const internal::kernel_retval< FullPivLUkernel () const
 
const internal::image_retval< FullPivLUimage (const MatrixType &originalMatrix) const
 
template<typename Rhs >
const internal::solve_retval< FullPivLU, Rhs > solve (const MatrixBase< Rhs > &b) 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 internal::solve_retval< FullPivLU, typename MatrixType::IdentityReturnType > inverse () const
 
MatrixType reconstructedMatrix () const
 
Index rows () const
 
Index cols () const
 

Protected Attributes

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

Detailed Description

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

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

Parameters
MatrixTypethe type of the matrix of which we are computing the LU decomposition

This class represents a LU decomposition of any matrix, with complete pivoting: the matrix A is decomposed as A = PLUQ 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:

Output:

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

Definition at line 216 of file ForwardDeclarations.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 387 of file FullPivLU.h.

388  : m_isInitialized(false), m_usePrescribedThreshold(false)
389 {
390 }
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 393 of file FullPivLU.h.

394  : m_lu(rows, cols),
395  m_p(rows),
396  m_q(cols),
397  m_rowsTranspositions(rows),
398  m_colsTranspositions(cols),
399  m_isInitialized(false),
400  m_usePrescribedThreshold(false)
401 {
402 }
template<typename MatrixType >
Eigen::FullPivLU< MatrixType >::FullPivLU ( const MatrixType &  matrix)

Constructor.

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

Definition at line 405 of file FullPivLU.h.

406  : m_lu(matrix.rows(), matrix.cols()),
407  m_p(matrix.rows()),
408  m_q(matrix.cols()),
409  m_rowsTranspositions(matrix.rows()),
410  m_colsTranspositions(matrix.cols()),
411  m_isInitialized(false),
412  m_usePrescribedThreshold(false)
413 {
414  compute(matrix);
415 }
Definition: math3d.h:219
FullPivLU & compute(const MatrixType &matrix)
Definition: FullPivLU.h:418

Member Function Documentation

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

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

419 {
420  // the permutations are stored as int indices, so just to be sure:
421  eigen_assert(matrix.rows()<=NumTraits<int>::highest() && matrix.cols()<=NumTraits<int>::highest());
422 
423  m_isInitialized = true;
424  m_lu = matrix;
425 
426  const Index size = matrix.diagonalSize();
427  const Index rows = matrix.rows();
428  const Index cols = matrix.cols();
429 
430  // will store the transpositions, before we accumulate them at the end.
431  // can't accumulate on-the-fly because that will be done in reverse order for the rows.
432  m_rowsTranspositions.resize(matrix.rows());
433  m_colsTranspositions.resize(matrix.cols());
434  Index number_of_transpositions = 0; // number of NONTRIVIAL transpositions, i.e. m_rowsTranspositions[i]!=i
435 
436  m_nonzero_pivots = size; // the generic case is that in which all pivots are nonzero (invertible case)
437  m_maxpivot = RealScalar(0);
438 
439  for(Index k = 0; k < size; ++k)
440  {
441  // First, we need to find the pivot.
442 
443  // biggest coefficient in the remaining bottom-right corner (starting at row k, col k)
444  Index row_of_biggest_in_corner, col_of_biggest_in_corner;
445  RealScalar biggest_in_corner;
446  biggest_in_corner = m_lu.bottomRightCorner(rows-k, cols-k)
447  .cwiseAbs()
448  .maxCoeff(&row_of_biggest_in_corner, &col_of_biggest_in_corner);
449  row_of_biggest_in_corner += k; // correct the values! since they were computed in the corner,
450  col_of_biggest_in_corner += k; // need to add k to them.
451 
452  if(biggest_in_corner==RealScalar(0))
453  {
454  // before exiting, make sure to initialize the still uninitialized transpositions
455  // in a sane state without destroying what we already have.
456  m_nonzero_pivots = k;
457  for(Index i = k; i < size; ++i)
458  {
459  m_rowsTranspositions.coeffRef(i) = i;
460  m_colsTranspositions.coeffRef(i) = i;
461  }
462  break;
463  }
464 
465  if(biggest_in_corner > m_maxpivot) m_maxpivot = biggest_in_corner;
466 
467  // Now that we've found the pivot, we need to apply the row/col swaps to
468  // bring it to the location (k,k).
469 
470  m_rowsTranspositions.coeffRef(k) = row_of_biggest_in_corner;
471  m_colsTranspositions.coeffRef(k) = col_of_biggest_in_corner;
472  if(k != row_of_biggest_in_corner) {
473  m_lu.row(k).swap(m_lu.row(row_of_biggest_in_corner));
474  ++number_of_transpositions;
475  }
476  if(k != col_of_biggest_in_corner) {
477  m_lu.col(k).swap(m_lu.col(col_of_biggest_in_corner));
478  ++number_of_transpositions;
479  }
480 
481  // Now that the pivot is at the right location, we update the remaining
482  // bottom-right corner by Gaussian elimination.
483 
484  if(k<rows-1)
485  m_lu.col(k).tail(rows-k-1) /= m_lu.coeff(k,k);
486  if(k<size-1)
487  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);
488  }
489 
490  // the main loop is over, we still have to accumulate the transpositions to find the
491  // permutations P and Q
492 
493  m_p.setIdentity(rows);
494  for(Index k = size-1; k >= 0; --k)
495  m_p.applyTranspositionOnTheRight(k, m_rowsTranspositions.coeff(k));
496 
497  m_q.setIdentity(cols);
498  for(Index k = 0; k < size; ++k)
499  m_q.applyTranspositionOnTheRight(k, m_colsTranspositions.coeff(k));
500 
501  m_det_pq = (number_of_transpositions%2) ? -1 : 1;
502  return *this;
503 }
Definition: math3d.h:219
Derived & applyTranspositionOnTheRight(Index i, Index j)
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 506 of file FullPivLU.h.

507 {
508  eigen_assert(m_isInitialized && "LU is not initialized.");
509  eigen_assert(m_lu.rows() == m_lu.cols() && "You can't take the determinant of a non-square matrix!");
510  return Scalar(m_det_pq) * Scalar(m_lu.diagonal().prod());
511 }
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 311 of file FullPivLU.h.

312  {
313  eigen_assert(m_isInitialized && "LU is not initialized.");
314  return cols() - rank();
315  }
Index rank() const
Definition: FullPivLU.h:294
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 kernel.
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:

Output:

See also
kernel()

Definition at line 187 of file FullPivLU.h.

188  {
189  eigen_assert(m_isInitialized && "LU is not initialized.");
190  return internal::image_retval<FullPivLU>(*this, originalMatrix);
191  }
template<typename _MatrixType>
const internal::solve_retval<FullPivLU,typename MatrixType::IdentityReturnType> Eigen::FullPivLU< _MatrixType >::inverse ( void  ) 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 362 of file FullPivLU.h.

363  {
364  eigen_assert(m_isInitialized && "LU is not initialized.");
365  eigen_assert(m_lu.rows() == m_lu.cols() && "You can't take the inverse of a non-square matrix!");
366  return internal::solve_retval<FullPivLU,typename MatrixType::IdentityReturnType>
367  (*this, MatrixType::Identity(m_lu.rows(), m_lu.cols()));
368  }
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 324 of file FullPivLU.h.

325  {
326  eigen_assert(m_isInitialized && "LU is not initialized.");
327  return rank() == cols();
328  }
Index rank() const
Definition: FullPivLU.h:294
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 349 of file FullPivLU.h.

350  {
351  eigen_assert(m_isInitialized && "LU is not initialized.");
352  return isInjective() && (m_lu.rows() == m_lu.cols());
353  }
bool isInjective() const
Definition: FullPivLU.h:324
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 337 of file FullPivLU.h.

338  {
339  eigen_assert(m_isInitialized && "LU is not initialized.");
340  return rank() == rows();
341  }
Index rank() const
Definition: FullPivLU.h:294
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:

Output:

See also
image()

Definition at line 161 of file FullPivLU.h.

162  {
163  eigen_assert(m_isInitialized && "LU is not initialized.");
164  return internal::kernel_retval<FullPivLU>(*this);
165  }
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 103 of file FullPivLU.h.

104  {
105  eigen_assert(m_isInitialized && "LU is not initialized.");
106  return m_lu;
107  }
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 125 of file FullPivLU.h.

125 { return m_maxpivot; }
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 116 of file FullPivLU.h.

117  {
118  eigen_assert(m_isInitialized && "LU is not initialized.");
119  return m_nonzero_pivots;
120  }
template<typename _MatrixType>
const PermutationPType& Eigen::FullPivLU< _MatrixType >::permutationP ( ) const
inline
Returns
the permutation matrix P
See also
permutationQ()

Definition at line 131 of file FullPivLU.h.

132  {
133  eigen_assert(m_isInitialized && "LU is not initialized.");
134  return m_p;
135  }
template<typename _MatrixType>
const PermutationQType& Eigen::FullPivLU< _MatrixType >::permutationQ ( ) const
inline
Returns
the permutation matrix Q
See also
permutationP()

Definition at line 141 of file FullPivLU.h.

142  {
143  eigen_assert(m_isInitialized && "LU is not initialized.");
144  return m_q;
145  }
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 294 of file FullPivLU.h.

295  {
296  using std::abs;
297  eigen_assert(m_isInitialized && "LU is not initialized.");
298  RealScalar premultiplied_threshold = abs(m_maxpivot) * threshold();
299  Index result = 0;
300  for(Index i = 0; i < m_nonzero_pivots; ++i)
301  result += (abs(m_lu.coeff(i,i)) > premultiplied_threshold);
302  return result;
303  }
RealScalar threshold() const
Definition: FullPivLU.h:279
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 purpose.

Definition at line 517 of file FullPivLU.h.

518 {
519  eigen_assert(m_isInitialized && "LU is not initialized.");
520  const Index smalldim = (std::min)(m_lu.rows(), m_lu.cols());
521  // LU
522  MatrixType res(m_lu.rows(),m_lu.cols());
523  // FIXME the .toDenseMatrix() should not be needed...
524  res = m_lu.leftCols(smalldim)
525  .template triangularView<UnitLower>().toDenseMatrix()
526  * m_lu.topRows(smalldim)
527  .template triangularView<Upper>().toDenseMatrix();
528 
529  // P^{-1}(LU)
530  res = m_p.inverse() * res;
531 
532  // (P^{-1}LU)Q^{-1}
533  res = res * m_q.inverse();
534 
535  return res;
536 }
Transpose< PermutationBase > inverse() const
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 254 of file FullPivLU.h.

255  {
256  m_usePrescribedThreshold = true;
257  m_prescribedThreshold = threshold;
258  return *this;
259  }
RealScalar threshold() const
Definition: FullPivLU.h:279
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 269 of file FullPivLU.h.

270  {
271  m_usePrescribedThreshold = false;
272  return *this;
273  }
template<typename _MatrixType>
template<typename Rhs >
const internal::solve_retval<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:

Output:

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

Definition at line 214 of file FullPivLU.h.

215  {
216  eigen_assert(m_isInitialized && "LU is not initialized.");
217  return internal::solve_retval<FullPivLU, Rhs>(*this, b.derived());
218  }
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 279 of file FullPivLU.h.

280  {
281  eigen_assert(m_isInitialized || m_usePrescribedThreshold);
282  return m_usePrescribedThreshold ? m_prescribedThreshold
283  // this formula comes from experimenting (see "LU precision tuning" thread on the list)
284  // and turns out to be identical to Higham's formula used already in LDLt.
285  : NumTraits<Scalar>::epsilon() * m_lu.diagonalSize();
286  }

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