Eigne Matrix Class Library
Dependents: Eigen_test Odometry_test AttitudeEstimation_usingTicker MPU9250_Quaternion_Binary_Serial ... more
FullPivLU< _MatrixType > Class Template Reference
[LU module]
LU decomposition of a matrix with complete pivoting, and related features. More...
#include <FullPivLU.h>
Public Member Functions | |
FullPivLU () | |
Default Constructor. | |
FullPivLU (Index rows, Index cols) | |
Default Constructor with memory preallocation. | |
FullPivLU (const MatrixType &matrix) | |
Constructor. | |
FullPivLU & | compute (const MatrixType &matrix) |
Computes the LU decomposition of the given matrix. | |
const MatrixType & | matrixLU () const |
Index | nonzeroPivots () const |
RealScalar | maxPivot () const |
const PermutationPType & | permutationP () const |
const PermutationQType & | permutationQ () const |
const internal::kernel_retval < FullPivLU > | kernel () const |
const internal::image_retval < FullPivLU > | image (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 |
FullPivLU & | setThreshold (const RealScalar &threshold) |
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. | |
FullPivLU & | setThreshold (Default_t) |
Allows to come back to the default behavior, letting Eigen use its default formula for determining the threshold. | |
RealScalar | threshold () const |
Returns the threshold that will be used by certain methods such as rank(). | |
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 |
Detailed Description
template<typename _MatrixType>
class Eigen::FullPivLU< _MatrixType >
LU decomposition of a matrix with complete pivoting, and related features.
- 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 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:
Definition at line 46 of file FullPivLU.h.
Constructor & Destructor Documentation
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 394 of file FullPivLU.h.
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 400 of file FullPivLU.h.
FullPivLU | ( | const MatrixType & | matrix ) |
Constructor.
- Parameters:
-
matrix the matrix of which to compute the LU decomposition. It is required to be nonzero.
Definition at line 412 of file FullPivLU.h.
Member Function Documentation
FullPivLU< MatrixType > & compute | ( | const MatrixType & | matrix ) |
Computes the LU decomposition of the given matrix.
- Parameters:
-
matrix the matrix of which to compute the LU decomposition. It is required to be nonzero.
- Returns:
- a reference to *this
Definition at line 425 of file FullPivLU.h.
internal::traits< MatrixType >::Scalar 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 515 of file FullPivLU.h.
Index dimensionOfKernel | ( | ) | const |
- 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 312 of file FullPivLU.h.
const internal::image_retval<FullPivLU> image | ( | const MatrixType & | originalMatrix ) | const |
- 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:
-
originalMatrix the 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 188 of file FullPivLU.h.
const internal::solve_retval<FullPivLU,typename MatrixType::IdentityReturnType> inverse | ( | ) | const |
- 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 363 of file FullPivLU.h.
bool isInjective | ( | ) | const |
- 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 325 of file FullPivLU.h.
bool isInvertible | ( | ) | const |
- 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 350 of file FullPivLU.h.
bool isSurjective | ( | ) | const |
- 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 338 of file FullPivLU.h.
const internal::kernel_retval<FullPivLU> kernel | ( | ) | const |
- 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 162 of file FullPivLU.h.
const MatrixType& matrixLU | ( | ) | const |
- 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 104 of file FullPivLU.h.
RealScalar maxPivot | ( | ) | const |
- Returns:
- the absolute value of the biggest pivot, i.e. the biggest diagonal coefficient of U.
Definition at line 126 of file FullPivLU.h.
Index nonzeroPivots | ( | ) | const |
- 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 117 of file FullPivLU.h.
const PermutationPType& permutationP | ( | ) | const |
- Returns:
- the permutation matrix P
- See also:
- permutationQ()
Definition at line 132 of file FullPivLU.h.
const PermutationQType& permutationQ | ( | ) | const |
- Returns:
- the permutation matrix Q
- See also:
- permutationP()
Definition at line 142 of file FullPivLU.h.
Index rank | ( | ) | const |
- 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 295 of file FullPivLU.h.
MatrixType reconstructedMatrix | ( | ) | const |
- Returns:
- the matrix represented by the decomposition, i.e., it returns the product: . This function is provided for debug purposes.
Definition at line 526 of file FullPivLU.h.
FullPivLU& setThreshold | ( | Default_t | ) |
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 270 of file FullPivLU.h.
FullPivLU& setThreshold | ( | const RealScalar & | threshold ) |
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:
-
threshold The new value to use as the threshold.
A pivot will be considered nonzero if its absolute value is strictly greater than where maxpivot is the biggest pivot.
If you want to come back to the default behavior, call setThreshold(Default_t)
Definition at line 255 of file FullPivLU.h.
const internal::solve_retval<FullPivLU, Rhs> solve | ( | const MatrixBase< Rhs > & | b ) | const |
- Returns:
- a solution x to the equation Ax=b, where A is the matrix of which *this is the LU decomposition.
- Parameters:
-
b the 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 215 of file FullPivLU.h.
RealScalar threshold | ( | ) | const |
Returns the threshold that will be used by certain methods such as rank().
See the documentation of setThreshold(const RealScalar&).
Definition at line 280 of file FullPivLU.h.
Generated on Tue Jul 12 2022 17:47:06 by 1.7.2