mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-01-24 14:45:14 +08:00
bug #899: make sparseqr unit test more stable by 1) trying with larger threshold and 2) relax rank computation for rank-deficient problems.
This commit is contained in:
parent
482c5fb321
commit
3b5deeb546
@ -43,6 +43,7 @@ int generate_sparse_rectangular_problem(MatrixType& A, DenseMat& dA, int maxRows
|
|||||||
|
|
||||||
template<typename Scalar> void test_sparseqr_scalar()
|
template<typename Scalar> void test_sparseqr_scalar()
|
||||||
{
|
{
|
||||||
|
typedef typename NumTraits<Scalar>::Real RealScalar;
|
||||||
typedef SparseMatrix<Scalar,ColMajor> MatrixType;
|
typedef SparseMatrix<Scalar,ColMajor> MatrixType;
|
||||||
typedef Matrix<Scalar,Dynamic,Dynamic> DenseMat;
|
typedef Matrix<Scalar,Dynamic,Dynamic> DenseMat;
|
||||||
typedef Matrix<Scalar,Dynamic,1> DenseVector;
|
typedef Matrix<Scalar,Dynamic,1> DenseVector;
|
||||||
@ -91,14 +92,34 @@ template<typename Scalar> void test_sparseqr_scalar()
|
|||||||
exit(0);
|
exit(0);
|
||||||
return;
|
return;
|
||||||
}
|
}
|
||||||
|
|
||||||
VERIFY_IS_APPROX(A * x, b);
|
// Compare with a dense QR solver
|
||||||
|
|
||||||
//Compare with a dense QR solver
|
|
||||||
ColPivHouseholderQR<DenseMat> dqr(dA);
|
ColPivHouseholderQR<DenseMat> dqr(dA);
|
||||||
refX = dqr.solve(b);
|
refX = dqr.solve(b);
|
||||||
|
|
||||||
VERIFY_IS_EQUAL(dqr.rank(), solver.rank());
|
bool rank_deficient = A.cols()>A.rows() || dqr.rank()<A.cols();
|
||||||
|
if(rank_deficient)
|
||||||
|
{
|
||||||
|
// rank deficient problem -> we might have to increase the threshold
|
||||||
|
// to get a correct solution.
|
||||||
|
RealScalar th = RealScalar(20)*dA.colwise().norm().maxCoeff()*(A.rows()+A.cols()) * NumTraits<RealScalar>::epsilon();
|
||||||
|
for(Index k=0; (k<16) && !test_isApprox(A*x,b); ++k)
|
||||||
|
{
|
||||||
|
th *= RealScalar(10);
|
||||||
|
solver.setPivotThreshold(th);
|
||||||
|
solver.compute(A);
|
||||||
|
x = solver.solve(b);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
VERIFY_IS_APPROX(A * x, b);
|
||||||
|
|
||||||
|
// For rank deficient problem, the estimated rank might
|
||||||
|
// be slightly off, so let's only raise a warning in such cases.
|
||||||
|
if(rank_deficient) ++g_test_level;
|
||||||
|
VERIFY_IS_EQUAL(solver.rank(), dqr.rank());
|
||||||
|
if(rank_deficient) --g_test_level;
|
||||||
|
|
||||||
if(solver.rank()==A.cols()) // full rank
|
if(solver.rank()==A.cols()) // full rank
|
||||||
VERIFY_IS_APPROX(x, refX);
|
VERIFY_IS_APPROX(x, refX);
|
||||||
// else
|
// else
|
||||||
|
Loading…
Reference in New Issue
Block a user