LU decomposition of a matrix with complete pivoting, and related features. More...
#include <FullPivLU.h>
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:
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; Eigen::FullPivLU<Matrix5x3> lu(m); 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;
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
typedef MatrixType::Index Index |
typedef internal::plain_col_type<MatrixType, Index>::type IntColVectorType |
typedef internal::plain_row_type<MatrixType, Index>::type IntRowVectorType |
typedef _MatrixType MatrixType |
typedef NumTraits<typename MatrixType::Scalar>::Real RealScalar |
typedef MatrixType::Scalar Scalar |
typedef internal::traits<MatrixType>::StorageKind StorageKind |
anonymous enum |
FullPivLU | ( | ) |
Default Constructor.
The default constructor is useful in cases in which the user intends to perform decompositions via LU::compute(const MatrixType&).
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 FullPivLU< _MatrixType >::compute().
References FullPivLU< _MatrixType >::m_lu.
Referenced by FullPivLU< _MatrixType >::dimensionOfKernel(), and FullPivLU< _MatrixType >::isInjective().
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< _MatrixType >::rows().
Referenced by FullPivLU< _MatrixType >::FullPivLU().
internal::traits< MatrixType >::Scalar determinant | ( | ) | const |
Index dimensionOfKernel | ( | ) | const [inline] |
References FullPivLU< _MatrixType >::cols(), FullPivLU< _MatrixType >::m_isInitialized, and FullPivLU< _MatrixType >::rank().
const internal::image_retval<FullPivLU> image | ( | const MatrixType & | originalMatrix | ) | const [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:
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:
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
References FullPivLU< _MatrixType >::m_isInitialized.
const internal::solve_retval<FullPivLU,typename MatrixType::IdentityReturnType> inverse | ( | ) | const [inline] |
References FullPivLU< _MatrixType >::m_isInitialized, and FullPivLU< _MatrixType >::m_lu.
bool isInjective | ( | ) | const [inline] |
References FullPivLU< _MatrixType >::cols(), FullPivLU< _MatrixType >::m_isInitialized, and FullPivLU< _MatrixType >::rank().
Referenced by FullPivLU< _MatrixType >::isInvertible().
bool isInvertible | ( | ) | const [inline] |
References FullPivLU< _MatrixType >::isInjective(), FullPivLU< _MatrixType >::m_isInitialized, and FullPivLU< _MatrixType >::m_lu.
bool isSurjective | ( | ) | const [inline] |
References FullPivLU< _MatrixType >::m_isInitialized, FullPivLU< _MatrixType >::rank(), and FullPivLU< _MatrixType >::rows().
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:
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: -9.46e-09 1.96e-08 2.53e-10 3.82e-08 6.06e-10 1.45e-08
References FullPivLU< _MatrixType >::m_isInitialized.
const MatrixType& matrixLU | ( | ) | const [inline] |
References FullPivLU< _MatrixType >::m_isInitialized, and FullPivLU< _MatrixType >::m_lu.
RealScalar maxPivot | ( | ) | const [inline] |
References FullPivLU< _MatrixType >::m_maxpivot.
Index nonzeroPivots | ( | ) | const [inline] |
References FullPivLU< _MatrixType >::m_isInitialized, and FullPivLU< _MatrixType >::m_nonzero_pivots.
const PermutationPType& permutationP | ( | ) | const [inline] |
References FullPivLU< _MatrixType >::m_isInitialized, and FullPivLU< _MatrixType >::m_p.
const PermutationQType& permutationQ | ( | ) | const [inline] |
References FullPivLU< _MatrixType >::m_isInitialized, and FullPivLU< _MatrixType >::m_q.
References abs(), FullPivLU< _MatrixType >::m_isInitialized, FullPivLU< _MatrixType >::m_lu, FullPivLU< _MatrixType >::m_maxpivot, FullPivLU< _MatrixType >::m_nonzero_pivots, and FullPivLU< _MatrixType >::threshold().
Referenced by FullPivLU< _MatrixType >::dimensionOfKernel(), FullPivLU< _MatrixType >::isInjective(), and FullPivLU< _MatrixType >::isSurjective().
MatrixType reconstructedMatrix | ( | ) | const |
References FullPivLU< _MatrixType >::m_lu.
Referenced by FullPivLU< _MatrixType >::compute(), and FullPivLU< _MatrixType >::isSurjective().
FullPivLU& 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.
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)
References FullPivLU< _MatrixType >::m_prescribedThreshold, FullPivLU< _MatrixType >::m_usePrescribedThreshold, and FullPivLU< _MatrixType >::threshold().
FullPivLU& 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&).
References FullPivLU< _MatrixType >::m_usePrescribedThreshold.
const internal::solve_retval<FullPivLU, Rhs> solve | ( | const MatrixBase< Rhs > & | b | ) | const [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:
bool a_solution_exists = (A*result).isApprox(b, 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(); Matrix2f y = Matrix2f::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;
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
References FullPivLU< _MatrixType >::m_isInitialized.
RealScalar threshold | ( | ) | const [inline] |
Returns the threshold that will be used by certain methods such as rank().
See the documentation of setThreshold(const RealScalar&).
References FullPivLU< _MatrixType >::m_isInitialized, FullPivLU< _MatrixType >::m_lu, FullPivLU< _MatrixType >::m_prescribedThreshold, and FullPivLU< _MatrixType >::m_usePrescribedThreshold.
Referenced by FullPivLU< _MatrixType >::rank(), and FullPivLU< _MatrixType >::setThreshold().
IntRowVectorType m_colsTranspositions [protected] |
bool m_isInitialized [protected] |
Referenced by FullPivLU< _MatrixType >::dimensionOfKernel(), FullPivLU< _MatrixType >::image(), FullPivLU< _MatrixType >::inverse(), FullPivLU< _MatrixType >::isInjective(), FullPivLU< _MatrixType >::isInvertible(), FullPivLU< _MatrixType >::isSurjective(), FullPivLU< _MatrixType >::kernel(), FullPivLU< _MatrixType >::matrixLU(), FullPivLU< _MatrixType >::nonzeroPivots(), FullPivLU< _MatrixType >::permutationP(), FullPivLU< _MatrixType >::permutationQ(), FullPivLU< _MatrixType >::rank(), FullPivLU< _MatrixType >::solve(), and FullPivLU< _MatrixType >::threshold().
MatrixType m_lu [protected] |
RealScalar m_maxpivot [protected] |
Referenced by FullPivLU< _MatrixType >::maxPivot(), and FullPivLU< _MatrixType >::rank().
Index m_nonzero_pivots [protected] |
Referenced by FullPivLU< _MatrixType >::nonzeroPivots(), and FullPivLU< _MatrixType >::rank().
PermutationPType m_p [protected] |
Referenced by FullPivLU< _MatrixType >::permutationP().
RealScalar m_prescribedThreshold [protected] |
Referenced by FullPivLU< _MatrixType >::setThreshold(), and FullPivLU< _MatrixType >::threshold().
PermutationQType m_q [protected] |
Referenced by FullPivLU< _MatrixType >::permutationQ().
IntColVectorType m_rowsTranspositions [protected] |
bool m_usePrescribedThreshold [protected] |
Referenced by FullPivLU< _MatrixType >::setThreshold(), and FullPivLU< _MatrixType >::threshold().