Symmlq.cxx
1 // Copyright (C) 2003-2009 Marc DuruflĂ©
2 //
3 // This file is part of the linear-algebra library Seldon,
4 // http://seldon.sourceforge.net/.
5 //
6 // Seldon is free software; you can redistribute it and/or modify it under the
7 // terms of the GNU Lesser General Public License as published by the Free
8 // Software Foundation; either version 2.1 of the License, or (at your option)
9 // any later version.
10 //
11 // Seldon is distributed in the hope that it will be useful, but WITHOUT ANY
12 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
13 // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
14 // more details.
15 //
16 // You should have received a copy of the GNU Lesser General Public License
17 // along with Seldon. If not, see http://www.gnu.org/licenses/.
18 
19 
20 #ifndef SELDON_FILE_ITERATIVE_SYMMLQ_CXX
21 
22 namespace Seldon
23 {
24 
26 
45 #ifdef SELDON_WITH_VIRTUAL
46  template<class T, class Vector1>
47  int Symmlq(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 Symmlq(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 alpha, beta, ibeta, beta_old, beta1,
62  ceta, ceta_oold, ceta_old, ceta_bar;
63  Complexe c, cold, s, sold, coold, soold, rho0, rho1, rho2, rho3, dp;
64  Complexe zero, one;
65  SetComplexZero(zero);
66  SetComplexOne(one);
67  ceta = zero;
68 
69  Vector1 r(b), z(b), u(b), v(b), w(b), u_old(b), v_old(b), w_bar(b);
70 
71  typedef typename ClassComplexType<Complexe>::Treal Treal;
72  Treal np, s_prod;
73  u_old.Fill(zero); v_old.Fill(zero); w.Fill(zero); w_bar.Fill(zero);
74 
75  int success_init = iter.Init(b);
76  if (success_init != 0)
77  return iter.ErrorCode();
78 
79  Copy(b, r);
80  // r = b - A x
81  if (!iter.IsInitGuess_Null())
82  iter.MltAdd(-one, A, x, one, r);
83  else
84  x.Fill(zero);
85 
86  ceta_oold = zero; ceta_old = zero;
87  c = one; cold = one; s = zero; sold = zero;
88  M.Solve(A, r, z);
89 
90  dp = DotProd(r, z);
91  dp = sqrt(dp); beta = dp; beta1 = beta;
92  s_prod = abs(beta1);
93 
94  Copy(r, v); Copy(z, u);
95  ibeta = one/beta;
96  Mlt(ibeta, v); Mlt(ibeta, u);
97  Copy(u, w_bar);
98  np = Norm2(b);
99 
100  iter.SetNumberIteration(0);
101  // Loop until the stopping criteria are satisfied
102  while (!iter.Finished(np))
103  {
104  // update
105  if (!iter.First())
106  {
107  Copy(v, v_old); Copy(u, u_old);
108  ibeta = one/beta;
109  // v = ibeta r
110  // u = ibeta z
111  Copy(r, v); Mlt(ibeta, v);
112  Copy(z, u); Mlt(ibeta, u);
113  // w = c*w_bar + s*u
114  Copy(w_bar, w); Mlt(c, w); Add(s, u, w);
115  // w_bar = -s*w_bar + c*u
116  Mlt(-s,w_bar); Add(c, u, w_bar);
117  // x = x+ceta*w
118  Add(ceta, w, x);
119 
120  ceta_oold = ceta_old;
121  ceta_old = ceta;
122  }
123 
124  // product matrix vector r = A u
125  iter.Mlt(A, u, r);
126  alpha = DotProd(u, r);
127  // preconditioning
128  M.Solve(A, r, z);
129 
130  // r = r - alpha*v
131  // z = z - alpha*u
132  Add(-alpha, v, r);
133  Add(-alpha, u, z);
134 
135  // r = r - beta*v_old
136  // z = z - beta*u_old
137  Add(-beta, v_old, r);
138  Add(-beta, u_old, z);
139 
140  beta_old = beta;
141  dp = DotProd(r, z);
142  beta = sqrt(dp);
143 
144  // QR factorization
145  coold = cold; cold = c; soold = sold; sold = s;
146  rho0 = cold * alpha - coold * sold * beta_old; // gamma_bar
147  rho1 = sqrt(rho0*rho0 + beta*beta); // gamma
148  rho2 = sold * alpha + coold * cold * beta_old; // delta
149  rho3 = soold * beta_old; // epsilon
150 
151  // Givens rotation
152  c = rho0 / rho1; s = beta / rho1;
153 
154  if (iter.First())
155  ceta = beta1/rho1;
156  else
157  ceta = -(rho2*ceta_old + rho3*ceta_oold)/rho1;
158 
159  s_prod *= abs(s);
160  if (c == zero)
161  np = s_prod*1e16;
162  else
163  np = s_prod/abs(c);
164 
165  ++iter;
166  }
167  if (c == zero)
168  ceta_bar = ceta*1e15;
169  else
170  ceta_bar = ceta/c;
171 
172  Add(ceta_bar, w_bar, x);
173 
174  return iter.ErrorCode();
175  }
176 
177 } // end namespace
178 
179 #define SELDON_FILE_ITERATIVE_SYMMLQ_CXX
180 #endif
Seldon::Iteration::First
bool First() const
Returns true if it is the first iteration.
Definition: Iterative.cxx:236
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::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::Symmlq
int Symmlq(const Matrix1 &A, Vector1 &x, const Vector1 &b, Preconditioner &M, Iteration< Titer > &iter)
Solves a linear system by using Symmetric LQ (SymmLQ)
Definition: Symmlq.cxx:52
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