Cgs.cxx
1 // Copyright (C) 2003-2009 Marc DuruflĂ©
2 // Copyright (C) 2001-2009 Vivien Mallet
3 //
4 // This file is part of the linear-algebra library Seldon,
5 // http://seldon.sourceforge.net/.
6 //
7 // Seldon is free software; you can redistribute it and/or modify it under the
8 // terms of the GNU Lesser General Public License as published by the Free
9 // Software Foundation; either version 2.1 of the License, or (at your option)
10 // any later version.
11 //
12 // Seldon is distributed in the hope that it will be useful, but WITHOUT ANY
13 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
14 // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
15 // more details.
16 //
17 // You should have received a copy of the GNU Lesser General Public License
18 // along with Seldon. If not, see http://www.gnu.org/licenses/.
19 
20 
21 #ifndef SELDON_FILE_ITERATIVE_CGS_CXX
22 
23 namespace Seldon
24 {
25 
27 
45 #ifdef SELDON_WITH_VIRTUAL
46  template<class T, class Vector1>
47  int Cgs(const VirtualMatrix<T>& A, Vector1& x, const Vector1& b,
48  Preconditioner_Base<T>& M,
49  Iteration<typename ClassComplexType<T>::Treal>& iter)
50 #else
51  template <class Titer, class Matrix1, class Vector1, class Preconditioner>
52  int Cgs(const Matrix1& A, Vector1& x, const Vector1& b,
53  Preconditioner& M, Iteration<Titer> & iter)
54 #endif
55  {
56  const int N = A.GetM();
57  if (N <= 0)
58  return 0;
59 
60  typedef typename Vector1::value_type Complexe;
61  Complexe rho_1, rho_2, alpha, beta, delta;
62  Vector1 p(b), phat(b), q(b), qhat(b), vhat(b), u(b), uhat(b),
63  r(b), rtilde(b);
64  Complexe zero, one;
65  SetComplexZero(zero);
66  SetComplexOne(one);
67  rho_2 = zero;
68 
69  // we initialize iter
70  int success_init = iter.Init(b);
71  if (success_init != 0)
72  return iter.ErrorCode();
73 
74  // we compute the initial residual r = b - Ax
75  Copy(b,r);
76  if (!iter.IsInitGuess_Null())
77  iter.MltAdd(-one, A, x, one, r);
78  else
79  x.Fill(zero);
80 
81  Copy(r, rtilde);
82 
83  iter.SetNumberIteration(0);
84  // Loop until the stopping criteria are reached
85  while (! iter.Finished(r))
86  {
87  rho_1 = DotProd(rtilde, r);
88 
89  if (rho_1 == zero)
90  {
91  iter.Fail(1, "Cgs breakdown #1");
92  break;
93  }
94 
95  if (iter.First())
96  {
97  Copy(r, u);
98  Copy(u, p);
99  }
100  else
101  {
102  // u = r + beta*q
103  // p = beta*(beta*p +q) + u where beta = rho_i/rho_{i-1}
104  beta = rho_1 / rho_2;
105  Copy(r, u);
106  Add(beta, q, u);
107  Mlt(beta, p);
108  Add(one, q, p);
109  Mlt(beta, p);
110  Add(one, u, p);
111  }
112 
113  // preconditioning phat = M^{-1} p
114  M.Solve(A, p, phat);
115 
116  // matrix vector product vhat = A*phat
117  iter.Mlt(A, phat, vhat); ++iter;
118  delta = DotProd(rtilde, vhat);
119  if (delta == zero)
120  {
121  iter.Fail(2, "Cgs breakdown #2");
122  break;
123  }
124  // q = u-alpha*vhat where alpha = rho_i/(rtilde,vhat)
125  alpha = rho_1 /delta;
126  Copy(u, q);
127  Add(-alpha, vhat, q);
128 
129  // u =u+q
130  Add(one, q, u);
131  M.Solve(A, u, uhat);
132 
133  Add(alpha, uhat, x);
134  iter.Mlt(A, uhat, qhat);
135  Add(-alpha, qhat, r);
136 
137  rho_2 = rho_1;
138 
139  ++iter;
140  }
141 
142  return iter.ErrorCode();
143  }
144 
145 } // end namespace
146 
147 #define SELDON_FILE_ITERATIVE_CGS_CXX
148 #endif
149 
Seldon::Iteration::First
bool First() const
Returns true if it is the first iteration.
Definition: Iterative.cxx:236
Seldon::Cgs
int Cgs(const Matrix1 &A, Vector1 &x, const Vector1 &b, Preconditioner &M, Iteration< Titer > &iter)
Solves linear system using Conjugate Gradient Squared (CGS)
Definition: Cgs.cxx:52
Seldon::Iteration::Init
int Init(const Vector1 &r)
Initialization with the right hand side.
Definition: Iterative.cxx:218
Seldon::Iteration::Mlt
void Mlt(const Matrix1 &A, const Vector1 &x, Vector1 &y)
Computes y = A x.
Definition: IterativeInline.cxx:123
Seldon::Iteration::ErrorCode
int ErrorCode()
Returns the error code (if an error occured)
Definition: Iterative.cxx:351
Seldon::Iteration::Fail
void Fail(int i, const string &s)
Informs of a failure in the iterative solver.
Definition: Iterative.cxx:330
Seldon::Iteration::MltAdd
void MltAdd(const T1 &alpha, const Matrix1 &A, const Vector1 &x, const T1 &beta, Vector1 &y)
Computes y = beta y + alpha A x.
Definition: IterativeInline.cxx:109
Seldon::Iteration::SetNumberIteration
void SetNumberIteration(int nb)
Changes the number of iterations.
Definition: Iterative.cxx:177
Seldon
Seldon namespace.
Definition: Array.cxx:24
Seldon::Iteration
Class containing parameters for an iterative resolution.
Definition: Iterative.hxx:59
Seldon::Iteration::IsInitGuess_Null
bool IsInitGuess_Null() const
Returns true if the initial guess is null.
Definition: Iterative.cxx:247
Seldon::Iteration::Finished
bool Finished(const Vector1 &r) const
Returns true if the iterative solver has reached its end.
Definition: Iterative.cxx:264