aboutsummaryrefslogtreecommitdiffhomepage
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2015-10-26 16:16:24 +0100
committerGravatar Gael Guennebaud <g.gael@free.fr>2015-10-26 16:16:24 +0100
commita5324a131f3816c8312e27a9dc928b8d56d8cf3b (patch)
tree9d241029b473b1c833b6dc2eaf8ef23b170b05ae
parentf93654ae16b261e462ee00c5255072f8dd7d387b (diff)
bug #1092: fix iterative solver ctors for expressions as input
-rw-r--r--Eigen/src/IterativeLinearSolvers/BiCGSTAB.h3
-rw-r--r--Eigen/src/IterativeLinearSolvers/ConjugateGradient.h3
-rw-r--r--Eigen/src/IterativeLinearSolvers/LeastSquareConjugateGradient.h3
-rw-r--r--test/sparse_solver.h65
-rw-r--r--unsupported/Eigen/src/IterativeSolvers/DGMRES.h4
-rw-r--r--unsupported/Eigen/src/IterativeSolvers/GMRES.h3
-rw-r--r--unsupported/Eigen/src/IterativeSolvers/MINRES.h3
7 files changed, 52 insertions, 32 deletions
diff --git a/Eigen/src/IterativeLinearSolvers/BiCGSTAB.h b/Eigen/src/IterativeLinearSolvers/BiCGSTAB.h
index 76e86a94a..4be00da47 100644
--- a/Eigen/src/IterativeLinearSolvers/BiCGSTAB.h
+++ b/Eigen/src/IterativeLinearSolvers/BiCGSTAB.h
@@ -182,7 +182,8 @@ public:
* this class becomes invalid. Call compute() to update it with the new
* 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() {}
diff --git a/Eigen/src/IterativeLinearSolvers/ConjugateGradient.h b/Eigen/src/IterativeLinearSolvers/ConjugateGradient.h
index 59092dc18..dbedf28fd 100644
--- a/Eigen/src/IterativeLinearSolvers/ConjugateGradient.h
+++ b/Eigen/src/IterativeLinearSolvers/ConjugateGradient.h
@@ -185,7 +185,8 @@ public:
* this class becomes invalid. Call compute() to update it with the new
* 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() {}
diff --git a/Eigen/src/IterativeLinearSolvers/LeastSquareConjugateGradient.h b/Eigen/src/IterativeLinearSolvers/LeastSquareConjugateGradient.h
index b578b2a7f..1593c57b5 100644
--- a/Eigen/src/IterativeLinearSolvers/LeastSquareConjugateGradient.h
+++ b/Eigen/src/IterativeLinearSolvers/LeastSquareConjugateGradient.h
@@ -175,7 +175,8 @@ public:
* this class becomes invalid. Call compute() to update it with the new
* 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() {}
diff --git a/test/sparse_solver.h b/test/sparse_solver.h
index a0254ff1c..b67653496 100644
--- a/test/sparse_solver.h
+++ b/test/sparse_solver.h
@@ -63,32 +63,47 @@ void check_sparse_solving(Solver& solver, const typename Solver::MatrixType& A,
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());
- Solver solver2(A);
- 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());
- DenseRhs oldb(db);
- x.setZero();
- 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;
- A2.reserve((ArrayXf::Random(A.outerSize())+2).template cast<typename Mat::StorageIndex>().eval());
- solver.compute(A2);
- Rhs x = solver.solve(b);
- VERIFY(x.isApprox(refX,test_precision<Scalar>()));
+ // test initialization ctor
+ {
+ Rhs x(b.rows(), b.cols());
+ Solver solver2(A);
+ 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());
+ DenseRhs oldb(db);
+ x.setZero();
+ 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;
+ A2.reserve((ArrayXf::Random(A.outerSize())+2).template cast<typename Mat::StorageIndex>().eval());
+ 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>()));
+ }
}
}
diff --git a/unsupported/Eigen/src/IterativeSolvers/DGMRES.h b/unsupported/Eigen/src/IterativeSolvers/DGMRES.h
index 52eb65a2f..ab82e782d 100644
--- a/unsupported/Eigen/src/IterativeSolvers/DGMRES.h
+++ b/unsupported/Eigen/src/IterativeSolvers/DGMRES.h
@@ -134,8 +134,8 @@ class DGMRES : public IterativeSolverBase<DGMRES<_MatrixType,_Preconditioner> >
* this class becomes invalid. Call compute() to update it with the new
* 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() {}
diff --git a/unsupported/Eigen/src/IterativeSolvers/GMRES.h b/unsupported/Eigen/src/IterativeSolvers/GMRES.h
index 05e5862a5..2cfa60140 100644
--- a/unsupported/Eigen/src/IterativeSolvers/GMRES.h
+++ b/unsupported/Eigen/src/IterativeSolvers/GMRES.h
@@ -288,7 +288,8 @@ public:
* this class becomes invalid. Call compute() to update it with the new
* 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() {}
diff --git a/unsupported/Eigen/src/IterativeSolvers/MINRES.h b/unsupported/Eigen/src/IterativeSolvers/MINRES.h
index c393112a4..84e491fa1 100644
--- a/unsupported/Eigen/src/IterativeSolvers/MINRES.h
+++ b/unsupported/Eigen/src/IterativeSolvers/MINRES.h
@@ -227,7 +227,8 @@ namespace Eigen {
* this class becomes invalid. Call compute() to update it with the new
* 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. */
~MINRES(){}