QCgs.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_QCGS_CXX
22 
23 namespace Seldon
24 {
25 
27 
46 #ifdef SELDON_WITH_VIRTUAL
47  template<class T, class Vector1>
48  int QCgs(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 QCgs(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  Complexe rho_1, rho_2, mu, nu, alpha, beta, sigma, delta;
63  Vector1 p(b), q(b), r(b), rtilde(b), u(b),
64  phat(b), r_qcgs(b), x_qcgs(b), v(b);
65 
66  Complexe zero, one;
67  SetComplexZero(zero);
68  SetComplexOne(one);
69 
70  // we initialize iter
71  int success_init = iter.Init(b);
72  if (success_init != 0)
73  return iter.ErrorCode();
74 
75  // we compute the residual r = b - Ax
76  Copy(b,r);
77  if (!iter.IsInitGuess_Null())
78  iter.MltAdd(-one, A, x, one, r);
79  else
80  x.Fill(zero);
81 
82  Copy(r, rtilde);
83  rho_1 = one;
84  q.Fill(zero); p.Fill(zero);
85  Copy(r, r_qcgs); Copy(x, x_qcgs);
86  Matrix<Complexe, Symmetric, RowSymPacked> bt_b(2,2), bt_b_m1(2,2);
87  Vector<Complexe> bt_rn(2);
88 
89  iter.SetNumberIteration(0);
90  // Loop until the stopping criteria are reached
91  while (! iter.Finished(r_qcgs))
92  {
93  rho_2 = DotProd(rtilde, r);
94 
95  if (rho_1 == zero)
96  {
97  iter.Fail(1, "QCgs breakdown #1");
98  break;
99  }
100  beta = rho_2/rho_1;
101 
102  // u = r + beta*q
103  // p = beta*(beta*p +q) + u where beta = rho_i/rho_{i-1}
104  Copy(r, u);
105  Add(beta, q, u);
106  Mlt(beta, p);
107  Add(one, q, p);
108  Mlt(beta, p);
109  Add(one, u, p);
110 
111  // preconditioning phat = M^{-1} p
112  M.Solve(A, p, phat);
113 
114  // product matrix vector vhat = A*phat
115  iter.Mlt(A, phat, v); ++iter;
116  sigma = DotProd(rtilde, v);
117  if (sigma == zero)
118  {
119  iter.Fail(2, "Qcgs breakdown #2");
120  break;
121  }
122  // q = u-alpha*v where alpha = rho_i/(rtilde,v)
123  alpha = rho_2 /sigma;
124  Copy(u, q);
125  Add(-alpha, v, q);
126 
127  // u = u +q
128  Add(one, q, u);
129  // x = x + alpha u
130  Add(alpha, u, x);
131  // preconditioning phat = M^{-1} u
132  M.Solve(A, u, phat);
133  // product matrix vector q = A*phat
134  iter.Mlt(A, phat, u);
135 
136  // r = r - alpha u
137  Add(-alpha, u, r);
138  // u = r_qcgs - r_n
139  Copy(r_qcgs, u);
140  Add(-one, r, u);
141  bt_b(0,0) = DotProd(u, u);
142  bt_b(1,1) = DotProd(v, v);
143  bt_b(1,0) = DotProd(u, v);
144 
145  // we compute inverse of bt_b
146  delta = bt_b(0,0)*bt_b(1,1) - bt_b(1,0)*bt_b(0,1);
147  if (delta == zero)
148  {
149  iter.Fail(3, "Qcgs breakdown #3");
150  break;
151  }
152  bt_b_m1(0,0) = bt_b(1,1)/delta;
153  bt_b_m1(1,1) = bt_b(0,0)/delta;
154  bt_b_m1(1,0) = -bt_b(1,0)/delta;
155 
156  bt_rn(0) = -DotProd(u, r); bt_rn(1) = -DotProd(v, r);
157  mu = bt_b_m1(0,0)*bt_rn(0)+bt_b_m1(0,1)*bt_rn(1);
158  nu = bt_b_m1(1,0)*bt_rn(0)+bt_b_m1(1,1)*bt_rn(1);
159 
160  // smoothing step
161  // r_qcgs = r + mu (r_qcgs - r) + nu v
162  Copy(r, r_qcgs); Add(mu, u, r_qcgs);
163  Add(nu, v, r_qcgs);
164  // x_qcgs = x + mu (x_qcgs - x) - nu p
165  Mlt(mu, x_qcgs);
166  Add(one-mu, x, x_qcgs);
167  Add(-nu, p, x_qcgs);
168 
169  rho_1 = rho_2;
170 
171  ++iter;
172  }
173 
174  M.Solve(A, x_qcgs, x);
175  return iter.ErrorCode();
176  }
177 
178 } // end namespace
179 
180 #define SELDON_FILE_ITERATIVE_QCGS_CXX
181 #endif
182 
Seldon::QCgs
int QCgs(const Matrix1 &A, Vector1 &x, const Vector1 &b, Preconditioner &M, Iteration< Titer > &iter)
Solves linear system using Quasi-minimized Conjugate Gradient Squared.
Definition: QCgs.cxx:53
Seldon::Vector
Definition: SeldonHeader.hxx:207
Seldon::Iteration::Init
int Init(const Vector1 &r)
Initialization with the right hand side.
Definition: Iterative.cxx:218
Seldon::Matrix
Definition: SeldonHeader.hxx:226
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