BiCg.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_BICG_CXX
22 
23 namespace Seldon
24 {
25 
27 
46 #ifdef SELDON_WITH_VIRTUAL
47  template<class T, class Vector1>
48  int BiCg(const VirtualMatrix<T>& A, Vector1& x, const Vector1& b,
49  Preconditioner_Base<T>& M,
50  Iteration<typename ClassComplexType<T>::Treal>& iter)
51 #else
52  template <class Titer, class Matrix1, class Vector1, class Preconditioner>
53  int BiCg(const Matrix1& A, Vector1& x, const Vector1& b,
54  Preconditioner& M, Iteration<Titer> & iter)
55 #endif
56  {
57  int N = A.GetM();
58  if (N <= 0)
59  return 0;
60 
61  typedef typename Vector1::value_type Complexe;
62  Complexe rho_1, rho_2, alpha, beta, delta;
63 
64  Vector1 r(b), z(b), p(b), q(b);
65  Vector1 r_tilde(b), z_tilde(b), p_tilde(b), q_tilde(b);
66  Complexe one, zero;
67  SetComplexZero(zero);
68  SetComplexOne(one);
69  rho_2 = zero;
70 
71  // we initialize iter
72  int success_init = iter.Init(b);
73  if (success_init != 0)
74  return iter.ErrorCode();
75 
76  // we compute the residual r = b - Ax
77  Copy(b, r);
78  if (!iter.IsInitGuess_Null())
79  iter.MltAdd(-one, A, x, one, r);
80  else
81  x.Fill(zero);
82 
83  Copy(r, r_tilde);
84 
85  iter.SetNumberIteration(0);
86  // Loop until the stopping criteria are reached
87  while (! iter.Finished(r))
88  {
89  // preconditioning z = M^{-1} r and z_tilde = M^{-t} r_tilde
90  M.Solve(A, r, z);
91  M.TransSolve(A, r_tilde, z_tilde);
92  // rho_1 = (z,r_tilde)
93  rho_1 = DotProd(z, r_tilde);
94 
95  if (rho_1 == zero)
96  {
97  iter.Fail(1, "Bicg breakdown #1");
98  break;
99  }
100 
101  if (iter.First())
102  {
103  Copy(z, p);
104  Copy(z_tilde, p_tilde);
105  }
106  else
107  {
108  // p=beta*p+z where beta = rho_i/rho_{i-1}
109  // p_tilde=beta*p_tilde+z_tilde
110  beta = rho_1 / rho_2;
111  Mlt(beta, p);
112  Add(one, z, p);
113  Mlt(beta, p_tilde);
114  Add(one, z_tilde, p_tilde);
115  }
116 
117  // we do the product matrix vector and transpose matrix vector
118  // q = A*p q_tilde = A^t p_tilde
119  iter.Mlt(A, p, q);
120  ++iter;
121  iter.Mlt(SeldonTrans, A, p_tilde, q_tilde);
122 
123  delta = DotProd(p_tilde, q);
124  if (delta == zero)
125  {
126  iter.Fail(2, "Bicg breakdown #2");
127  break;
128  }
129 
130  alpha = rho_1 / delta;
131 
132  // the new iterate x=x+alpha*p and residual r=r-alpha*q
133  // where alpha = rho_i/delta
134  Add(alpha, p, x);
135  Add(-alpha, q, r);
136  Add(-alpha, q_tilde, r_tilde);
137 
138  rho_2 = rho_1;
139 
140  ++iter;
141  }
142 
143  return iter.ErrorCode();
144  }
145 
146 }
147 
148 #define SELDON_FILE_ITERATIVE_BICG_CXX
149 #endif
Seldon::Iteration::First
bool First() const
Returns true if it is the first iteration.
Definition: Iterative.cxx:236
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::BiCg
int BiCg(const Matrix1 &A, Vector1 &x, const Vector1 &b, Preconditioner &M, Iteration< Titer > &iter)
Solves a linear system by using BiConjugate Gradient (BICG)
Definition: BiCg.cxx:53
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