LU decomposition of a matrix with complete pivoting, and related features. More...
#include <FullPivLU.h>
Public Member Functions | |
| FullPivLU & | compute (const MatrixType &matrix) |
| internal::traits< MatrixType >::Scalar | determinant () const |
| Index | dimensionOfKernel () const |
| FullPivLU () | |
| Default Constructor. | |
| FullPivLU (const MatrixType &matrix) | |
| FullPivLU (Index rows, Index cols) | |
| Default Constructor with memory preallocation. | |
| const internal::image_retval< FullPivLU > | image (const MatrixType &originalMatrix) const |
| const internal::solve_retval< FullPivLU, typename MatrixType::IdentityReturnType > | inverse () const |
| bool | isInjective () const |
| bool | isInvertible () const |
| bool | isSurjective () const |
| const internal::kernel_retval< FullPivLU > | kernel () const |
| const MatrixType & | matrixLU () const |
| RealScalar | maxPivot () const |
| Index | nonzeroPivots () const |
| const PermutationPType & | permutationP () const |
| const PermutationQType & | permutationQ () const |
| Index | rank () const |
| MatrixType | reconstructedMatrix () const |
| FullPivLU & | setThreshold (const RealScalar &threshold) |
| FullPivLU & | setThreshold (Default_t) |
| template<typename Rhs> | |
| const internal::solve_retval< FullPivLU, Rhs > | solve (const MatrixBase< Rhs > &b) const |
| RealScalar | threshold () const |
LU decomposition of a matrix with complete pivoting, and related features.
| 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 = 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:
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
| FullPivLU | ( | ) |
Default Constructor.
The default constructor is useful in cases in which the user intends to perform decompositions via LU::compute(const MatrixType&).
Referenced by compute(), setThreshold(), and setThreshold().
| 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.
| FullPivLU | ( | const MatrixType & | matrix | ) |
Constructor.
| matrix | the matrix of which to compute the LU decomposition. It is required to be nonzero. |
References compute().
| FullPivLU< MatrixType > & compute | ( | const MatrixType & | matrix | ) |
Computes the LU decomposition of the given matrix.
| matrix | the matrix of which to compute the LU decomposition. It is required to be nonzero. |
References FullPivLU().
Referenced by FullPivLU().
| internal::traits< MatrixType >::Scalar determinant | ( | ) | const |
|
inline |
|
inline |
| 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. |
Example:
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
|
inline |
|
inline |
References rank().
Referenced by isInvertible().
|
inline |
References isInjective().
|
inline |
References rank().
|
inline |
Example:
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:
-1.12e-08 1.49e-08
-1.4e-09 -4.05e-08
1.49e-08 -2.98e-08
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
References threshold().
Referenced by isInjective(), and isSurjective().
| MatrixType reconstructedMatrix | ( | ) | const |
|
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.
| threshold | The new value to use as the threshold. |
A pivot will be considered nonzero if its absolute value is strictly greater than 
If you want to come back to the default behavior, call setThreshold(Default_t)
References FullPivLU(), and threshold().
|
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.
See the documentation of setThreshold(const RealScalar&).
References FullPivLU().
|
inline |
| 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. |
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:
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:
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
|
inline |
Returns the threshold that will be used by certain methods such as rank().
See the documentation of setThreshold(const RealScalar&).
Referenced by rank(), and setThreshold().