Householder rank-revealing QR decomposition of a matrix with column-pivoting.
_MatrixType | the type of the matrix of which we are computing the QR decomposition |
This class performs a rank-revealing QR decomposition of a matrix A into matrices P, Q and R such that
\[ \mathbf{A} \, \mathbf{P} = \mathbf{Q} \, \mathbf{R} \]
by using Householder transformations. Here, P is a permutation matrix, Q a unitary matrix and R an upper triangular matrix.
This decomposition performs column pivoting in order to be rank-revealing and improve numerical stability. It is slower than HouseholderQR, and faster than FullPivHouseholderQR.
This class supports the inplace decomposition mechanism.
MatrixType::RealScalar | absDeterminant () const |
ColPivHouseholderQR () | |
Default Constructor. More... |
|
template<typename InputType > | |
ColPivHouseholderQR (const EigenBase< InputType > &matrix) | |
Constructs a QR factorization from a given matrix. More... |
|
template<typename InputType > | |
ColPivHouseholderQR (EigenBase< InputType > &matrix) | |
Constructs a QR factorization from a given matrix. More... |
|
ColPivHouseholderQR (Index rows, Index cols) | |
Default Constructor with memory preallocation. More... |
|
const PermutationType & | colsPermutation () const |
template<typename InputType > | |
ColPivHouseholderQR< MatrixType > & | compute (const EigenBase< InputType > &matrix) |
Index | dimensionOfKernel () const |
const HCoeffsType & | hCoeffs () const |
HouseholderSequenceType | householderQ () const |
ComputationInfo | info () const |
Reports whether the QR factorization was successful. More... |
|
const Inverse< ColPivHouseholderQR > | inverse () const |
bool | isInjective () const |
bool | isInvertible () const |
bool | isSurjective () const |
MatrixType::RealScalar | logAbsDeterminant () const |
const MatrixType & | matrixQR () const |
const MatrixType & | matrixR () const |
RealScalar | maxPivot () const |
Index | nonzeroPivots () const |
Index | rank () const |
ColPivHouseholderQR & | setThreshold (const RealScalar &threshold) |
ColPivHouseholderQR & | setThreshold (Default_t) |
template<typename Rhs > | |
const Solve< ColPivHouseholderQR, Rhs > | solve (const MatrixBase< Rhs > &b) const |
RealScalar | threshold () const |
Public Member Functions inherited from Eigen::SolverBase< ColPivHouseholderQR< _MatrixType > > | |
AdjointReturnType | adjoint () const |
ColPivHouseholderQR< _MatrixType > & | derived () |
const ColPivHouseholderQR< _MatrixType > & | derived () const |
const Solve< ColPivHouseholderQR< _MatrixType >, Rhs > | solve (const MatrixBase< Rhs > &b) const |
SolverBase () | |
ConstTransposeReturnType | transpose () const |
Public Member Functions inherited from Eigen::EigenBase< Derived > | |
EIGEN_CONSTEXPR Index | cols () const EIGEN_NOEXCEPT |
Derived & | derived () |
const Derived & | derived () const |
EIGEN_CONSTEXPR Index | rows () const EIGEN_NOEXCEPT |
EIGEN_CONSTEXPR Index | size () const EIGEN_NOEXCEPT |
Public Types inherited from Eigen::EigenBase< Derived > | |
typedef Eigen::Index | Index |
The interface type of indices. More... |
|
| inline |
Default Constructor.
The default constructor is useful in cases in which the user intends to perform decompositions via ColPivHouseholderQR::compute(const MatrixType&).
| inline |
Default Constructor with memory preallocation.
Like the default constructor but with preallocation of the internal data according to the specified problem size.
| inlineexplicit |
Constructs a QR factorization from a given matrix.
This constructor computes the QR factorization of the matrix matrix by calling the method compute(). It is a short cut for:
ColPivHouseholderQR<MatrixType> qr(matrix.rows(), matrix.cols()); qr.compute(matrix);
| inlineexplicit |
Constructs a QR factorization from a given matrix.
This overloaded constructor is provided for inplace decomposition when MatrixType
is a Eigen::Ref.
MatrixType::RealScalar Eigen::ColPivHouseholderQR< MatrixType >::absDeterminant |
| inline |
ColPivHouseholderQR<MatrixType>& Eigen::ColPivHouseholderQR< _MatrixType >::compute | ( | const EigenBase< InputType > & | matrix | ) |
Performs the QR factorization of the given matrix matrix. The result of the factorization is stored into *this
, and a reference to *this
is returned.
| inline |
| inline |
Q
.For advanced uses only.
ColPivHouseholderQR< MatrixType >::HouseholderSequenceType Eigen::ColPivHouseholderQR< MatrixType >::householderQ |
qr.householderQ().setLength(qr.nonzeroPivots())
| inline |
Reports whether the QR factorization was successful.
Success
. It is provided for compatibility with other factorization routines. Success
| inline |
| inline |
| inline |
| inline |
MatrixType::RealScalar Eigen::ColPivHouseholderQR< MatrixType >::logAbsDeterminant |
| inline |
| inline |
| inline |
| inline |
| inline |
| 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 QR 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 \( \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)
| 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.
qr.setThreshold(Eigen::Default);
See the documentation of setThreshold(const RealScalar&).
| inline |
This method finds a solution x to the equation Ax=b, where A is the matrix of which *this is the QR decomposition, if any exists.
b | the right-hand-side of the equation to solve. |
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.
Example:
Matrix3f m = Matrix3f::Random(); Matrix3f y = Matrix3f::Random(); cout << "Here is the matrix m:" << endl << m << endl; cout << "Here is the matrix y:" << endl << y << endl; Matrix3f x; x = m.colPivHouseholderQr().solve(y); assert(y.isApprox(m*x)); cout << "Here is a solution x to the equation mx=y:" << endl << x << endl;
Output:
Here is the matrix m: 0.68 0.597 -0.33 -0.211 0.823 0.536 0.566 -0.605 -0.444 Here is the matrix y: 0.108 -0.27 0.832 -0.0452 0.0268 0.271 0.258 0.904 0.435 Here is a solution x to the equation mx=y: 0.609 2.68 1.67 -0.231 -1.57 0.0713 0.51 3.51 1.05
| inline |
Returns the threshold that will be used by certain methods such as rank().
See the documentation of setThreshold(const RealScalar&).
© Eigen.
Licensed under the MPL2 License.
https://eigen.tuxfamily.org/dox/classEigen_1_1ColPivHouseholderQR.html