MinRes.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_MINRES_CXX
22 
23 namespace Seldon
24 {
25 
27 
46 #ifdef SELDON_WITH_VIRTUAL
47  template<class T, class Vector1>
48  int MinRes(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 MinRes(const Matrix1& A, Vector1& x, const Vector1& b,
54  Preconditioner& M, Iteration<Titer> & iter)
55 #endif
56  {
57  const int N = A.GetM();
58  if (N <= 0)
59  return 0;
60 
61  typedef typename Vector1::value_type Complexe;
62  Vector1 u_old(b), u(b), r(b), v_old(b), v(b),
63  w_old(b), w(b), z(b), w_oold(b);
64 
65  Complexe dp, beta, ibeta, beta_old, alpha, eta, ceta;
66  Complexe cold, coold, c, soold, sold, s, rho0, rho1, rho2, rho3;
67  Complexe zero, one;
68  SetComplexZero(zero);
69  SetComplexOne(one);
70 
71  int success_init = iter.Init(b);
72  if (success_init != 0)
73  return iter.ErrorCode();
74 
75  Copy(b,r);
76  // r = b - A x
77  if (!iter.IsInitGuess_Null())
78  iter.MltAdd(-one, A, x, one, r);
79  else
80  x.Fill(zero);
81 
82  u_old.Fill(zero); v_old.Fill(zero);
83  w_old.Fill(zero); w.Fill(zero); w_oold.Fill(zero);
84 
85  // preconditioning
86  M.Solve(A, r, z);
87  dp = DotProd(r, z);
88  dp = sqrt(dp); beta = dp; eta = beta;
89  Copy(r, v); Copy(z, u);
90 
91  ibeta = one / beta;
92  Mlt(ibeta, v); Mlt(ibeta, u);
93 
94  c = one; s = zero; cold = one; sold = zero;
95  typedef typename ClassComplexType<Complexe>::Treal Treal;
96  Treal np = Norm2(b);
97 
98  iter.SetNumberIteration(0);
99  // Loop until the stopping criteria are satisfied
100  while (!iter.Finished(np))
101  {
102  // matrix-vector product r = A*u
103  iter.Mlt(A, u, r);
104  alpha = DotProd(r, u);
105  // preconditioning
106  M.Solve(A, r, z);
107 
108  // r = r - alpha v
109  // z = z - alpha u
110  Add(-alpha, v, r);
111  Add(-alpha, u, z);
112  // r = r - beta v_old
113  // z = z - beta u_old
114  Add(-beta, v_old, r);
115  Add(-beta, u_old, z);
116 
117  beta_old = beta;
118 
119  dp = DotProd(r, z);
120  beta = sqrt(dp);
121 
122  // QR factorization
123  coold = cold; cold = c; soold = sold; sold = s;
124 
125  rho0 = cold * alpha - coold * sold * beta_old;
126  rho1 = sqrt(rho0*rho0 + beta*beta);
127  rho2 = sold * alpha + coold * cold * beta_old;
128  rho3 = soold * beta_old;
129 
130  // Givens rotation
131  if (rho1 == zero)
132  {
133  iter.Fail(1, "Minres breakdown #1");
134  break;
135  }
136  c = rho0 / rho1;
137  s = beta / rho1;
138 
139  // update
140  Copy(w_old, w_oold); Copy(w, w_old);
141  Copy(u, w);
142 
143  Add(-rho2, w_old, w);
144  Add(-rho3, w_oold, w);
145  Mlt(one/rho1, w);
146 
147  ceta = c*eta;
148  Add(ceta, w, x);
149  eta = -s*eta;
150 
151  Copy(v, v_old); Copy(u, u_old);
152  Copy(r, v); Copy(z, u);
153  if (beta == zero)
154  {
155  iter.Fail(2, "MinRes breakdown #2");
156  break;
157  }
158  ibeta = one/beta;
159  Mlt(ibeta, v); Mlt(ibeta, u);
160 
161  // residual norm
162  np *= abs(s);
163  ++ iter;
164  }
165  return iter.ErrorCode();
166  }
167 
168 } // end namespace
169 
170 #define SELDON_FILE_ITERATIVE_MINRES_CXX
171 #endif
Seldon::Norm2
ClassComplexType< T >::Treal Norm2(const VectorExpression< T, E > &X)
returns 2-norm of an expression with vectors
Definition: Functions_BaseInline.cxx:153
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::MinRes
int MinRes(const Matrix1 &A, Vector1 &x, const Vector1 &b, Preconditioner &M, Iteration< Titer > &iter)
Solves a linear system by using Minimum Residual (MinRes)
Definition: MinRes.cxx:53
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