diff --git a/Eigen/src/IterativeLinearSolvers/IterativeSolverBase.h b/Eigen/src/IterativeLinearSolvers/IterativeSolverBase.h index 38334028a..46bc0ac78 100644 --- a/Eigen/src/IterativeLinearSolvers/IterativeSolverBase.h +++ b/Eigen/src/IterativeLinearSolvers/IterativeSolverBase.h @@ -54,9 +54,10 @@ public: */ template explicit IterativeSolverBase(const SparseMatrixBase& A) + : mp_matrix(A) { init(); - compute(A.derived()); + compute(mp_matrix); } ~IterativeSolverBase() {} @@ -69,7 +70,7 @@ public: template Derived& analyzePattern(const SparseMatrixBase& A) { - grab(A); + grab(A.derived()); m_preconditioner.analyzePattern(mp_matrix); m_isInitialized = true; m_analysisIsOk = true; @@ -90,7 +91,7 @@ public: Derived& factorize(const SparseMatrixBase& A) { eigen_assert(m_analysisIsOk && "You must first call analyzePattern()"); - grab(A); + grab(A.derived()); m_preconditioner.factorize(mp_matrix); m_factorizationIsOk = true; m_info = Success; @@ -110,7 +111,7 @@ public: template Derived& compute(const SparseMatrixBase& A) { - grab(A); + grab(A.derived()); m_preconditioner.compute(mp_matrix); m_isInitialized = true; m_analysisIsOk = true; @@ -229,6 +230,15 @@ protected: ::new (&mp_matrix) Ref(A); } + void grab(const Ref &A) + { + if(&(A.derived()) != &mp_matrix) + { + mp_matrix.~Ref(); + ::new (&mp_matrix) Ref(A); + } + } + MatrixType m_dummy; Ref mp_matrix; Preconditioner m_preconditioner; diff --git a/test/sparse_solver.h b/test/sparse_solver.h index 42b365eaa..45cfdad25 100644 --- a/test/sparse_solver.h +++ b/test/sparse_solver.h @@ -83,6 +83,15 @@ void check_sparse_solving(Solver& solver, const typename Solver::MatrixType& A, VERIFY(xm.isApprox(refX,test_precision())); } + // 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())); + } + // test dense Block as the result and rhs: { DenseRhs x(db.rows(), db.cols());