aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/QR
diff options
context:
space:
mode:
authorGravatar Benoit Jacob <jacob.benoit.1@gmail.com>2009-08-24 11:11:41 -0400
committerGravatar Benoit Jacob <jacob.benoit.1@gmail.com>2009-08-24 11:11:41 -0400
commit0eb142f5595aa7d18b6c08a9e8ebc355f3a9b525 (patch)
tree7f7ae202d86074cc8a93e7ff266f2f9f21cfa87f /Eigen/src/QR
parent3288e5157a8d2c8a35c5c0835e4670386cded0ff (diff)
bring the modern comfort also to ColPivotingHouseholderQR
+ some fixes in FullPivotingHouseholderQR
Diffstat (limited to 'Eigen/src/QR')
-rw-r--r--Eigen/src/QR/ColPivotingHouseholderQR.h159
-rw-r--r--Eigen/src/QR/FullPivotingHouseholderQR.h27
2 files changed, 171 insertions, 15 deletions
diff --git a/Eigen/src/QR/ColPivotingHouseholderQR.h b/Eigen/src/QR/ColPivotingHouseholderQR.h
index ed4b84f63..0aec6a607 100644
--- a/Eigen/src/QR/ColPivotingHouseholderQR.h
+++ b/Eigen/src/QR/ColPivotingHouseholderQR.h
@@ -31,14 +31,14 @@
*
* \class ColPivotingHouseholderQR
*
- * \brief Householder rank-revealing QR decomposition of a matrix
+ * \brief Householder rank-revealing QR decomposition of a matrix with column-pivoting
*
* \param MatrixType the type of the matrix of which we are computing the QR decomposition
*
* This class performs a rank-revealing QR decomposition using Householder transformations.
*
- * This decomposition performs full-pivoting in order to be rank-revealing and achieve optimal
- * numerical stability.
+ * This decomposition performs column pivoting in order to be rank-revealing and improve
+ * numerical stability. It is slower than HouseholderQR, and faster than FullPivotingHouseholderQR.
*
* \sa MatrixBase::colPivotingHouseholderQr()
*/
@@ -82,6 +82,8 @@ template<typename MatrixType> class ColPivotingHouseholderQR
/** 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.
*
+ * \returns \c true if a solution exists, \c false if no solution exists.
+ *
* \param b the right-hand-side of the equation to solve.
*
* \param result a pointer to the vector/matrix in which to store the solution, if any exists.
@@ -95,7 +97,7 @@ template<typename MatrixType> class ColPivotingHouseholderQR
* Output: \verbinclude ColPivotingHouseholderQR_solve.out
*/
template<typename OtherDerived, typename ResultType>
- void solve(const MatrixBase<OtherDerived>& b, ResultType *result) const;
+ bool solve(const MatrixBase<OtherDerived>& b, ResultType *result) const;
MatrixType matrixQ(void) const;
@@ -111,12 +113,122 @@ template<typename MatrixType> class ColPivotingHouseholderQR
return m_cols_permutation;
}
+ /** \returns the absolute value of the determinant of the matrix of which
+ * *this is the QR decomposition. It has only linear complexity
+ * (that is, O(n) where n is the dimension of the square matrix)
+ * as the QR decomposition has already been computed.
+ *
+ * \note This is only for square matrices.
+ *
+ * \warning a determinant can be very big or small, so for matrices
+ * of large enough dimension, there is a risk of overflow/underflow.
+ * One way to work around that is to use logAbsDeterminant() instead.
+ *
+ * \sa logAbsDeterminant(), MatrixBase::determinant()
+ */
+ typename MatrixType::RealScalar absDeterminant() const;
+
+ /** \returns the natural log of the absolute value of the determinant of the matrix of which
+ * *this is the QR decomposition. It has only linear complexity
+ * (that is, O(n) where n is the dimension of the square matrix)
+ * as the QR decomposition has already been computed.
+ *
+ * \note This is only for square matrices.
+ *
+ * \note This method is useful to work around the risk of overflow/underflow that's inherent
+ * to determinant computation.
+ *
+ * \sa absDeterminant(), MatrixBase::determinant()
+ */
+ typename MatrixType::RealScalar logAbsDeterminant() const;
+
+ /** \returns the rank of the matrix of which *this is the QR decomposition.
+ *
+ * \note This is computed at the time of the construction of the QR decomposition. This
+ * method does not perform any further computation.
+ */
inline int rank() const
{
ei_assert(m_isInitialized && "ColPivotingHouseholderQR is not initialized.");
return m_rank;
}
+ /** \returns the dimension of the kernel of the matrix of which *this is the QR decomposition.
+ *
+ * \note Since the rank is computed at the time of the construction of the QR decomposition, this
+ * method almost does not perform any further computation.
+ */
+ inline int dimensionOfKernel() const
+ {
+ ei_assert(m_isInitialized && "ColPivotingHouseholderQR is not initialized.");
+ return m_qr.cols() - m_rank;
+ }
+
+ /** \returns true if the matrix of which *this is the QR decomposition represents an injective
+ * linear map, i.e. has trivial kernel; false otherwise.
+ *
+ * \note Since the rank is computed at the time of the construction of the QR decomposition, this
+ * method almost does not perform any further computation.
+ */
+ inline bool isInjective() const
+ {
+ ei_assert(m_isInitialized && "ColPivotingHouseholderQR is not initialized.");
+ return m_rank == m_qr.cols();
+ }
+
+ /** \returns true if the matrix of which *this is the QR decomposition represents a surjective
+ * linear map; false otherwise.
+ *
+ * \note Since the rank is computed at the time of the construction of the QR decomposition, this
+ * method almost does not perform any further computation.
+ */
+ inline bool isSurjective() const
+ {
+ ei_assert(m_isInitialized && "ColPivotingHouseholderQR is not initialized.");
+ return m_rank == m_qr.rows();
+ }
+
+ /** \returns true if the matrix of which *this is the QR decomposition is invertible.
+ *
+ * \note Since the rank is computed at the time of the construction of the QR decomposition, this
+ * method almost does not perform any further computation.
+ */
+ inline bool isInvertible() const
+ {
+ ei_assert(m_isInitialized && "ColPivotingHouseholderQR is not initialized.");
+ return isInjective() && isSurjective();
+ }
+
+ /** Computes the inverse of the matrix of which *this is the QR decomposition.
+ *
+ * \param result a pointer to the matrix into which to store the inverse. Resized if needed.
+ *
+ * \note If this matrix is not invertible, *result is left with undefined coefficients.
+ * Use isInvertible() to first determine whether this matrix is invertible.
+ *
+ * \sa inverse()
+ */
+ inline void computeInverse(MatrixType *result) const
+ {
+ ei_assert(m_isInitialized && "ColPivotingHouseholderQR is not initialized.");
+ ei_assert(m_qr.rows() == m_qr.cols() && "You can't take the inverse of a non-square matrix!");
+ solve(MatrixType::Identity(m_qr.rows(), m_qr.cols()), result);
+ }
+
+ /** \returns the inverse of the matrix of which *this is the QR decomposition.
+ *
+ * \note If this matrix is not invertible, the returned matrix has undefined coefficients.
+ * Use isInvertible() to first determine whether this matrix is invertible.
+ *
+ * \sa computeInverse()
+ */
+ inline MatrixType inverse() const
+ {
+ MatrixType result;
+ computeInverse(&result);
+ return result;
+ }
+
protected:
MatrixType m_qr;
HCoeffsType m_hCoeffs;
@@ -130,6 +242,22 @@ template<typename MatrixType> class ColPivotingHouseholderQR
#ifndef EIGEN_HIDE_HEAVY_CODE
template<typename MatrixType>
+typename MatrixType::RealScalar ColPivotingHouseholderQR<MatrixType>::absDeterminant() const
+{
+ ei_assert(m_isInitialized && "ColPivotingHouseholderQR is not initialized.");
+ ei_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
+ return ei_abs(m_qr.diagonal().prod());
+}
+
+template<typename MatrixType>
+typename MatrixType::RealScalar ColPivotingHouseholderQR<MatrixType>::logAbsDeterminant() const
+{
+ ei_assert(m_isInitialized && "ColPivotingHouseholderQR is not initialized.");
+ ei_assert(m_qr.rows() == m_qr.cols() && "You can't take the determinant of a non-square matrix!");
+ return m_qr.diagonal().cwise().abs().cwise().log().sum();
+}
+
+template<typename MatrixType>
ColPivotingHouseholderQR<MatrixType>& ColPivotingHouseholderQR<MatrixType>::compute(const MatrixType& matrix)
{
int rows = matrix.rows();
@@ -199,12 +327,23 @@ ColPivotingHouseholderQR<MatrixType>& ColPivotingHouseholderQR<MatrixType>::comp
template<typename MatrixType>
template<typename OtherDerived, typename ResultType>
-void ColPivotingHouseholderQR<MatrixType>::solve(
+bool ColPivotingHouseholderQR<MatrixType>::solve(
const MatrixBase<OtherDerived>& b,
ResultType *result
) const
{
ei_assert(m_isInitialized && "ColPivotingHouseholderQR is not initialized.");
+ result->resize(m_qr.cols(), b.cols());
+ if(m_rank==0)
+ {
+ if(b.squaredNorm() == RealScalar(0))
+ {
+ result->setZero();
+ return true;
+ }
+ else return false;
+ }
+
const int rows = m_qr.rows();
const int cols = b.cols();
ei_assert(b.rows() == rows);
@@ -219,13 +358,21 @@ void ColPivotingHouseholderQR<MatrixType>::solve(
.applyHouseholderOnTheLeft(m_qr.col(k).end(remainingSize-1), m_hCoeffs.coeff(k), &temp.coeffRef(0));
}
+ if(!isSurjective())
+ {
+ // is c is in the image of R ?
+ RealScalar biggest_in_upper_part_of_c = c.corner(TopLeft, m_rank, c.cols()).cwise().abs().maxCoeff();
+ RealScalar biggest_in_lower_part_of_c = c.corner(BottomLeft, rows-m_rank, c.cols()).cwise().abs().maxCoeff();
+ if(!ei_isMuchSmallerThan(biggest_in_lower_part_of_c, biggest_in_upper_part_of_c, m_precision))
+ return false;
+ }
m_qr.corner(TopLeft, m_rank, m_rank)
.template triangularView<UpperTriangular>()
.solveInPlace(c.corner(TopLeft, m_rank, c.cols()));
- result->resize(m_qr.cols(), b.cols());
for(int i = 0; i < m_rank; ++i) result->row(m_cols_permutation.coeff(i)) = c.row(i);
for(int i = m_rank; i < m_qr.cols(); ++i) result->row(m_cols_permutation.coeff(i)).setZero();
+ return true;
}
/** \returns the matrix Q */
diff --git a/Eigen/src/QR/FullPivotingHouseholderQR.h b/Eigen/src/QR/FullPivotingHouseholderQR.h
index 77a0abedc..77b664f6e 100644
--- a/Eigen/src/QR/FullPivotingHouseholderQR.h
+++ b/Eigen/src/QR/FullPivotingHouseholderQR.h
@@ -31,7 +31,7 @@
*
* \class FullPivotingHouseholderQR
*
- * \brief Householder rank-revealing QR decomposition of a matrix
+ * \brief Householder rank-revealing QR decomposition of a matrix with full pivoting
*
* \param MatrixType the type of the matrix of which we are computing the QR decomposition
*
@@ -62,12 +62,11 @@ template<typename MatrixType> class FullPivotingHouseholderQR
typedef Matrix<Scalar, 1, ColsAtCompileTime> RowVectorType;
typedef Matrix<Scalar, RowsAtCompileTime, 1> ColVectorType;
- /**
- * \brief Default Constructor.
- *
- * The default constructor is useful in cases in which the user intends to
- * perform decompositions via FullPivotingHouseholderQR::compute(const MatrixType&).
- */
+ /** \brief Default Constructor.
+ *
+ * The default constructor is useful in cases in which the user intends to
+ * perform decompositions via FullPivotingHouseholderQR::compute(const MatrixType&).
+ */
FullPivotingHouseholderQR() : m_qr(), m_hCoeffs(), m_isInitialized(false) {}
FullPivotingHouseholderQR(const MatrixType& matrix)
@@ -81,6 +80,8 @@ template<typename MatrixType> class FullPivotingHouseholderQR
/** 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.
*
+ * \returns \c true if a solution exists, \c false if no solution exists.
+ *
* \param b the right-hand-side of the equation to solve.
*
* \param result a pointer to the vector/matrix in which to store the solution, if any exists.
@@ -345,7 +346,16 @@ bool FullPivotingHouseholderQR<MatrixType>::solve(
) const
{
ei_assert(m_isInitialized && "FullPivotingHouseholderQR is not initialized.");
- if(m_rank==0) return false;
+ result->resize(m_qr.cols(), b.cols());
+ if(m_rank==0)
+ {
+ if(b.squaredNorm() == RealScalar(0))
+ {
+ result->setZero();
+ return true;
+ }
+ else return false;
+ }
const int rows = m_qr.rows();
const int cols = b.cols();
@@ -374,7 +384,6 @@ bool FullPivotingHouseholderQR<MatrixType>::solve(
.template triangularView<UpperTriangular>()
.solveInPlace(c.corner(TopLeft, m_rank, c.cols()));
- result->resize(m_qr.cols(), b.cols());
for(int i = 0; i < m_rank; ++i) result->row(m_cols_permutation.coeff(i)) = c.row(i);
for(int i = m_rank; i < m_qr.cols(); ++i) result->row(m_cols_permutation.coeff(i)).setZero();
return true;