mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-08-03 10:40:39 +08:00
bug #1092: fix iterative solver ctors for expressions as input
This commit is contained in:
parent
9055400f3d
commit
d0980c7706
@ -186,7 +186,8 @@ public:
|
||||
* this class becomes invalid. Call compute() to update it with the new
|
||||
* matrix A, or modify a copy of A.
|
||||
*/
|
||||
BiCGSTAB(const MatrixType& A) : Base(A) {}
|
||||
template<typename MatrixDerived>
|
||||
explicit BiCGSTAB(const EigenBase<MatrixDerived>& A) : Base(A.derived()) {}
|
||||
|
||||
~BiCGSTAB() {}
|
||||
|
||||
|
@ -176,7 +176,8 @@ public:
|
||||
* this class becomes invalid. Call compute() to update it with the new
|
||||
* matrix A, or modify a copy of A.
|
||||
*/
|
||||
ConjugateGradient(const MatrixType& A) : Base(A) {}
|
||||
template<typename MatrixDerived>
|
||||
explicit ConjugateGradient(const EigenBase<MatrixDerived>& A) : Base(A.derived()) {}
|
||||
|
||||
~ConjugateGradient() {}
|
||||
|
||||
|
@ -67,6 +67,22 @@ void check_sparse_solving(Solver& solver, const typename Solver::MatrixType& A,
|
||||
VERIFY(oldb.isApprox(db) && "sparse solver testing: the rhs should not be modified!");
|
||||
VERIFY(x.isApprox(refX,test_precision<Scalar>()));
|
||||
}
|
||||
|
||||
// if not too large, do some extra check:
|
||||
if(A.rows()<2000)
|
||||
{
|
||||
|
||||
// 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>()));
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
template<typename Solver, typename Rhs>
|
||||
|
@ -133,8 +133,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() {}
|
||||
|
||||
|
@ -285,7 +285,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() {}
|
||||
|
||||
|
@ -228,7 +228,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(){}
|
||||
|
Loading…
x
Reference in New Issue
Block a user