aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/SparseQR
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2019-02-19 22:52:15 +0100
committerGravatar Gael Guennebaud <g.gael@free.fr>2019-02-19 22:52:15 +0100
commit482c5fb321695f7992d3bb718b7f64f2feaf61d5 (patch)
tree426d0c13f470ca29803c066d5bc66cce8f19306f /Eigen/src/SparseQR
parent9ac1634fdff94bf18b534066eb0e3029ac182fe2 (diff)
bug #899: remove "rank-revealing" qualifier for SparseQR and warn that it is not always rank-revealing.
Diffstat (limited to 'Eigen/src/SparseQR')
-rw-r--r--Eigen/src/SparseQR/SparseQR.h22
1 files changed, 17 insertions, 5 deletions
diff --git a/Eigen/src/SparseQR/SparseQR.h b/Eigen/src/SparseQR/SparseQR.h
index 1a28389e8..d1fb96f5c 100644
--- a/Eigen/src/SparseQR/SparseQR.h
+++ b/Eigen/src/SparseQR/SparseQR.h
@@ -41,15 +41,16 @@ namespace internal {
/**
* \ingroup SparseQR_Module
* \class SparseQR
- * \brief Sparse left-looking rank-revealing QR factorization
+ * \brief Sparse left-looking QR factorization with numerical column pivoting
*
- * This class implements a left-looking rank-revealing QR decomposition
- * of sparse matrices. When a column has a norm less than a given tolerance
+ * This class implements a left-looking QR decomposition of sparse matrices
+ * with numerical column pivoting.
+ * When a column has a norm less than a given tolerance
* it is implicitly permuted to the end. The QR factorization thus obtained is
* given by A*P = Q*R where R is upper triangular or trapezoidal.
*
* P is the column permutation which is the product of the fill-reducing and the
- * rank-revealing permutations. Use colsPermutation() to get it.
+ * numerical permutations. Use colsPermutation() to get it.
*
* Q is the orthogonal matrix represented as products of Householder reflectors.
* Use matrixQ() to get an expression and matrixQ().adjoint() to get the adjoint.
@@ -64,6 +65,17 @@ namespace internal {
*
* \implsparsesolverconcept
*
+ * The numerical pivoting strategy and default threshold are the same as in SuiteSparse QR, and
+ * detailed in the following paper:
+ * <i>
+ * Tim Davis, "Algorithm 915, SuiteSparseQR: Multifrontal Multithreaded Rank-Revealing
+ * Sparse QR Factorization, ACM Trans. on Math. Soft. 38(1), 2011.
+ * </i>
+ * Even though it is qualified as "rank-revealing", this strategy might fail for some
+ * rank deficient problems. When this class is used to solve linear or least-square problems
+ * it is thus strongly recommended to check the accuracy of the computed solution. If it
+ * failed, it usually helps to increase the threshold with setPivotThreshold.
+ *
* \warning The input sparse matrix A must be in compressed mode (see SparseMatrix::makeCompressed()).
* \warning For complex matrices matrixQ().transpose() will actually return the adjoint matrix.
*
@@ -331,7 +343,7 @@ void SparseQR<MatrixType,OrderingType>::analyzePattern(const MatrixType& mat)
m_R.resize(m, n);
m_Q.resize(m, diagSize);
- // Allocate space for nonzero elements : rough estimation
+ // Allocate space for nonzero elements: rough estimation
m_R.reserve(2*mat.nonZeros()); //FIXME Get a more accurate estimation through symbolic factorization with the etree
m_Q.reserve(2*mat.nonZeros());
m_hcoeffs.resize(diagSize);