aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/IterativeLinearSolvers
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2015-02-16 17:05:10 +0100
committerGravatar Gael Guennebaud <g.gael@free.fr>2015-02-16 17:05:10 +0100
commit0f464d9d876d9bf6bdfe64d7fe47cd8e35ad3c2c (patch)
tree7b6b7ab03c3b51afe17c4abe6941abbb9e89fc05 /Eigen/src/IterativeLinearSolvers
parent470d26d5803caaf6e41dda5d3e864d8757e32f2d (diff)
bug #897: fix regression in BiCGSTAB(mat) ctor (an all other iterative solvers).
Add respective regression unit test.
Diffstat (limited to 'Eigen/src/IterativeLinearSolvers')
-rw-r--r--Eigen/src/IterativeLinearSolvers/IterativeSolverBase.h18
1 files changed, 14 insertions, 4 deletions
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<typename SparseMatrixDerived>
explicit IterativeSolverBase(const SparseMatrixBase<SparseMatrixDerived>& A)
+ : mp_matrix(A)
{
init();
- compute(A.derived());
+ compute(mp_matrix);
}
~IterativeSolverBase() {}
@@ -69,7 +70,7 @@ public:
template<typename SparseMatrixDerived>
Derived& analyzePattern(const SparseMatrixBase<SparseMatrixDerived>& 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<SparseMatrixDerived>& 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<typename SparseMatrixDerived>
Derived& compute(const SparseMatrixBase<SparseMatrixDerived>& 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<const MatrixType>(A);
}
+ void grab(const Ref<const MatrixType> &A)
+ {
+ if(&(A.derived()) != &mp_matrix)
+ {
+ mp_matrix.~Ref<const MatrixType>();
+ ::new (&mp_matrix) Ref<const MatrixType>(A);
+ }
+ }
+
MatrixType m_dummy;
Ref<const MatrixType> mp_matrix;
Preconditioner m_preconditioner;