21 #ifndef SELDON_FILE_ITERATIVE_CG_CXX
45 #ifdef SELDON_WITH_VIRTUAL
46 template<
class T,
class Vector1>
47 int Cg(
const VirtualMatrix<T>& A, Vector1& x,
const Vector1& b,
48 Preconditioner_Base<T>& M,
49 Iteration<
typename ClassComplexType<T>::Treal>& iter)
51 template <
class Titer,
class Matrix1,
class Vector1,
class Preconditioner>
52 int Cg(
const Matrix1& A, Vector1& x,
const Vector1& b,
56 const int N = A.GetM();
60 typedef typename Vector1::value_type Complexe;
61 Complexe rho, rho_1, alpha, beta, delta;
62 Vector1 p(b), q(b), r(b), z(b);
66 rho = one; rho_1 = one;
69 int success_init = iter.
Init(b);
70 if (success_init != 0)
76 iter.
MltAdd(-one, A, x, one, r);
89 rho = DotProdConj(r, z);
93 iter.
Fail(1,
"Cg breakdown #1");
109 delta = DotProdConj(p, q);
112 iter.
Fail(2,
"Cg breakdown #2");
132 #define SELDON_FILE_ITERATIVE_CG_CXX
bool First() const
Returns true if it is the first iteration.
int Init(const Vector1 &r)
Initialization with the right hand side.
void Mlt(const Matrix1 &A, const Vector1 &x, Vector1 &y)
Computes y = A x.
int Cg(const Matrix1 &A, Vector1 &x, const Vector1 &b, Preconditioner &M, Iteration< Titer > &iter)
Solves a linear system by using Conjugate Gradient (CG)
int ErrorCode()
Returns the error code (if an error occured)
void Fail(int i, const string &s)
Informs of a failure in the iterative solver.
void MltAdd(const T1 &alpha, const Matrix1 &A, const Vector1 &x, const T1 &beta, Vector1 &y)
Computes y = beta y + alpha A x.
void SetNumberIteration(int nb)
Changes the number of iterations.
Class containing parameters for an iterative resolution.
bool IsInitGuess_Null() const
Returns true if the initial guess is null.
bool Finished(const Vector1 &r) const
Returns true if the iterative solver has reached its end.