mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-05-01 08:14:10 +08:00
bug #941: fix accuracy issue in ColPivHouseholderQR, do not stop decomposition on a small pivot
This commit is contained in:
parent
9d82f7e30d
commit
f1092d2f73
@ -476,20 +476,10 @@ ColPivHouseholderQR<MatrixType>& ColPivHouseholderQR<MatrixType>::compute(const
|
||||
// we store that back into our table: it can't hurt to correct our table.
|
||||
m_colSqNorms.coeffRef(biggest_col_index) = biggest_col_sq_norm;
|
||||
|
||||
// if the current biggest column is smaller than epsilon times the initial biggest column,
|
||||
// terminate to avoid generating nan/inf values.
|
||||
// Note that here, if we test instead for "biggest == 0", we get a failure every 1000 (or so)
|
||||
// repetitions of the unit test, with the result of solve() filled with large values of the order
|
||||
// of 1/(size*epsilon).
|
||||
if(biggest_col_sq_norm < threshold_helper * RealScalar(rows-k))
|
||||
{
|
||||
// Track the number of meaningful pivots but do not stop the decomposition to make
|
||||
// sure that the initial matrix is properly reproduced. See bug 941.
|
||||
if(m_nonzero_pivots==size && biggest_col_sq_norm < threshold_helper * RealScalar(rows-k))
|
||||
m_nonzero_pivots = k;
|
||||
m_hCoeffs.tail(size-k).setZero();
|
||||
m_qr.bottomRightCorner(rows-k,cols-k)
|
||||
.template triangularView<StrictlyLower>()
|
||||
.setZero();
|
||||
break;
|
||||
}
|
||||
|
||||
// apply the transposition to the columns
|
||||
m_colsTranspositions.coeffRef(k) = biggest_col_index;
|
||||
@ -518,7 +508,7 @@ ColPivHouseholderQR<MatrixType>& ColPivHouseholderQR<MatrixType>::compute(const
|
||||
}
|
||||
|
||||
m_colsPermutation.setIdentity(PermIndexType(cols));
|
||||
for(PermIndexType k = 0; k < m_nonzero_pivots; ++k)
|
||||
for(PermIndexType k = 0; k < size/*m_nonzero_pivots*/; ++k)
|
||||
m_colsPermutation.applyTranspositionOnTheRight(k, PermIndexType(m_colsTranspositions.coeff(k)));
|
||||
|
||||
m_det_pq = (number_of_transpositions%2) ? -1 : 1;
|
||||
@ -574,13 +564,15 @@ struct Assignment<DstXprType, Inverse<ColPivHouseholderQR<MatrixType> >, interna
|
||||
|
||||
} // end namespace internal
|
||||
|
||||
/** \returns the matrix Q as a sequence of householder transformations */
|
||||
/** \returns the matrix Q as a sequence of householder transformations.
|
||||
* You can extract the meaningful part only by using:
|
||||
* \code qr.householderQ().setLength(qr.nonzeroPivots()) */
|
||||
template<typename MatrixType>
|
||||
typename ColPivHouseholderQR<MatrixType>::HouseholderSequenceType ColPivHouseholderQR<MatrixType>
|
||||
::householderQ() const
|
||||
{
|
||||
eigen_assert(m_isInitialized && "ColPivHouseholderQR is not initialized.");
|
||||
return HouseholderSequenceType(m_qr, m_hCoeffs.conjugate()).setLength(m_nonzero_pivots);
|
||||
return HouseholderSequenceType(m_qr, m_hCoeffs.conjugate());
|
||||
}
|
||||
|
||||
#ifndef __CUDACC__
|
||||
|
Loading…
x
Reference in New Issue
Block a user