bug #1092: fix iterative solver ctors for expressions as input

This commit is contained in:
Gael Guennebaud 2015-10-26 16:16:24 +01:00
parent f93654ae16
commit a5324a131f
7 changed files with 52 additions and 32 deletions

View File

@ -182,7 +182,8 @@ public:
* this class becomes invalid. Call compute() to update it with the new * this class becomes invalid. Call compute() to update it with the new
* matrix A, or modify a copy of A. * matrix A, or modify a copy of A.
*/ */
explicit BiCGSTAB(const MatrixType& A) : Base(A) {} template<typename MatrixDerived>
explicit BiCGSTAB(const EigenBase<MatrixDerived>& A) : Base(A.derived()) {}
~BiCGSTAB() {} ~BiCGSTAB() {}

View File

@ -185,7 +185,8 @@ public:
* this class becomes invalid. Call compute() to update it with the new * this class becomes invalid. Call compute() to update it with the new
* matrix A, or modify a copy of A. * matrix A, or modify a copy of A.
*/ */
explicit ConjugateGradient(const MatrixType& A) : Base(A) {} template<typename MatrixDerived>
explicit ConjugateGradient(const EigenBase<MatrixDerived>& A) : Base(A.derived()) {}
~ConjugateGradient() {} ~ConjugateGradient() {}

View File

@ -175,7 +175,8 @@ public:
* this class becomes invalid. Call compute() to update it with the new * this class becomes invalid. Call compute() to update it with the new
* matrix A, or modify a copy of A. * matrix A, or modify a copy of A.
*/ */
explicit LeastSquaresConjugateGradient(const MatrixType& A) : Base(A) {} template<typename MatrixDerived>
explicit LeastSquaresConjugateGradient(const EigenBase<MatrixDerived>& A) : Base(A.derived()) {}
~LeastSquaresConjugateGradient() {} ~LeastSquaresConjugateGradient() {}

View File

@ -63,32 +63,47 @@ void check_sparse_solving(Solver& solver, const typename Solver::MatrixType& A,
VERIFY(xm.isApprox(refX,test_precision<Scalar>())); VERIFY(xm.isApprox(refX,test_precision<Scalar>()));
} }
// test initialization ctor // if not too large, do some extra check:
if(A.rows()<2000)
{ {
Rhs x(b.rows(), b.cols()); // test initialization ctor
Solver solver2(A); {
VERIFY(solver2.info() == Success); Rhs x(b.rows(), b.cols());
x = solver2.solve(b); Solver solver2(A);
VERIFY(x.isApprox(refX,test_precision<Scalar>())); VERIFY(solver2.info() == Success);
} x = solver2.solve(b);
VERIFY(x.isApprox(refX,test_precision<Scalar>()));
// test dense Block as the result and rhs: }
{
DenseRhs x(refX.rows(), refX.cols()); // test dense Block as the result and rhs:
DenseRhs oldb(db); {
x.setZero(); DenseRhs x(refX.rows(), refX.cols());
x.block(0,0,x.rows(),x.cols()) = solver.solve(db.block(0,0,db.rows(),db.cols())); DenseRhs oldb(db);
VERIFY(oldb.isApprox(db) && "sparse solver testing: the rhs should not be modified!"); x.setZero();
VERIFY(x.isApprox(refX,test_precision<Scalar>())); x.block(0,0,x.rows(),x.cols()) = solver.solve(db.block(0,0,db.rows(),db.cols()));
} VERIFY(oldb.isApprox(db) && "sparse solver testing: the rhs should not be modified!");
VERIFY(x.isApprox(refX,test_precision<Scalar>()));
// test uncompressed inputs }
{
Mat A2 = A; // test uncompressed inputs
A2.reserve((ArrayXf::Random(A.outerSize())+2).template cast<typename Mat::StorageIndex>().eval()); {
solver.compute(A2); Mat A2 = A;
Rhs x = solver.solve(b); A2.reserve((ArrayXf::Random(A.outerSize())+2).template cast<typename Mat::StorageIndex>().eval());
VERIFY(x.isApprox(refX,test_precision<Scalar>())); solver.compute(A2);
Rhs x = solver.solve(b);
VERIFY(x.isApprox(refX,test_precision<Scalar>()));
}
// test expression as input
{
solver.compute(0.5*(A+A));
Rhs x = solver.solve(b);
VERIFY(x.isApprox(refX,test_precision<Scalar>()));
Solver solver2(0.5*(A+A));
Rhs x2 = solver2.solve(b);
VERIFY(x2.isApprox(refX,test_precision<Scalar>()));
}
} }
} }

View File

@ -134,8 +134,8 @@ class DGMRES : public IterativeSolverBase<DGMRES<_MatrixType,_Preconditioner> >
* this class becomes invalid. Call compute() to update it with the new * this class becomes invalid. Call compute() to update it with the new
* matrix A, or modify a copy of A. * matrix A, or modify a copy of A.
*/ */
DGMRES(const MatrixType& A) : Base(A),m_restart(30),m_neig(0),m_r(0),m_maxNeig(5),m_isDeflAllocated(false),m_isDeflInitialized(false) template<typename MatrixDerived>
{} explicit DGMRES(const EigenBase<MatrixDerived>& A) : Base(A.derived()), m_restart(30),m_neig(0),m_r(0),m_maxNeig(5),m_isDeflAllocated(false),m_isDeflInitialized(false) {}
~DGMRES() {} ~DGMRES() {}

View File

@ -288,7 +288,8 @@ public:
* this class becomes invalid. Call compute() to update it with the new * this class becomes invalid. Call compute() to update it with the new
* matrix A, or modify a copy of A. * matrix A, or modify a copy of A.
*/ */
GMRES(const MatrixType& A) : Base(A), m_restart(30) {} template<typename MatrixDerived>
explicit GMRES(const EigenBase<MatrixDerived>& A) : Base(A.derived()), m_restart(30) {}
~GMRES() {} ~GMRES() {}

View File

@ -227,7 +227,8 @@ namespace Eigen {
* this class becomes invalid. Call compute() to update it with the new * this class becomes invalid. Call compute() to update it with the new
* matrix A, or modify a copy of A. * matrix A, or modify a copy of A.
*/ */
MINRES(const MatrixType& A) : Base(A) {} template<typename MatrixDerived>
explicit MINRES(const EigenBase<MatrixDerived>& A) : Base(A.derived()) {}
/** Destructor. */ /** Destructor. */
~MINRES(){} ~MINRES(){}