OOPS
FGMRES.h
Go to the documentation of this file.
1 /*
2  * (C) Copyright 2009-2016 ECMWF.
3  *
4  * This software is licensed under the terms of the Apache Licence Version 2.0
5  * which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
6  * In applying this licence, ECMWF does not waive the privileges and immunities
7  * granted to it by virtue of its status as an intergovernmental organisation nor
8  * does it submit to any jurisdiction.
9  */
10 
11 #ifndef OOPS_ASSIMILATION_FGMRES_H_
12 #define OOPS_ASSIMILATION_FGMRES_H_
13 
14 #include <cmath>
15 #include <iostream>
16 #include <vector>
17 
21 #include "oops/util/dot_product.h"
22 #include "oops/util/formats.h"
23 #include "oops/util/Logger.h"
24 
25 namespace oops {
26 
27 /*! \file FGMRES.h
28  * \brief FGMRES solver for Ax=b.
29  *
30  * Flexible GMRES solver for Ax=b.(based on implementation following
31  * Youcef Saad,SIAM Journal on Scientific Computing, Volume 14 Issue 2,
32  * March 1993)
33 
34  * A must be square, but need not be symmetric.
35  * A preconditioner must be supplied that, given a
36  * vector q, returns an approximate solution of Ap=q.
37  *
38  * On entry:
39  * - x = starting point, \f$ X_0 \f$.
40  * - b = right hand side.
41  * - A = \f$ A \f$.
42  * - precond = preconditioner \f$ F_k \approx (A)^{-1} \f$.
43  * - maxiter = maximum number of iterations
44  * - tol = error tolerance
45  *
46  * On exit, x will contain the solution.
47  * The return value is the achieved reduction in residual norm.
48  *
49  * Iteration will stop if the maximum iteration limit "maxiter" is reached
50  * or if the residual norm reduces by a factor of "tolerance".
51  *
52  * VECTOR must implement:
53  * - dot_product
54  * - operator(=)
55  * - operator(+=),
56  * - operator(-=)
57  * - operator(*=) [double * VECTOR],
58  * - axpy
59  *
60  * AMATRIX and PMATRIX must implement a method:
61  * - void multiply(const VECTOR&, VECTOR&) const
62  *
63  * which applies the matrix to the first argument, and returns the
64  * matrix-vector product in the second. (Note: the const is optional, but
65  * recommended.)
66  */
67 
68 template <typename VECTOR, typename AMATRIX, typename PMATRIX>
69 double FGMRES(VECTOR & x, const VECTOR & b,
70  const AMATRIX & A, const PMATRIX & precond,
71  const int maxiter, const double tolerance) {
72  std::vector<VECTOR> V;
73  std::vector<VECTOR> Z;
74  std::vector< std::vector<double> > H;
75  std::vector<double> cs(maxiter+1, 0);
76  std::vector<double> sn(maxiter+1, 0);
77  std::vector<double> s;
78  std::vector<double> y(maxiter+1, 0);
79 
80  VECTOR r(x);
81  VECTOR z(x);
82  VECTOR work(x);
83  VECTOR work2(x);
84  VECTOR w(x);
85  VECTOR xinit(x);
86 
87  r = b;
88  double xnrm2 = dot_product(x, x);
89  if (xnrm2 > 0.0) {
90  A.multiply(x, work);
91  r -= work; // r = b - Ax
92  }
93 
94  double bnrm2 = sqrt(dot_product(b, b));
95  double rnrm2 = sqrt(dot_product(r, r));
96  double normReduction = rnrm2 / bnrm2;
97 
98  // Initialiaze (maxiter + 1) by maxiter matrix H
99  H.resize(maxiter);
100  for (int ii = 0; ii <= maxiter-1; ii++) {
101  H[ii].resize(maxiter + 1);
102  for (int jj = 0; jj <= maxiter; jj++) {
103  H[ii][jj] = 0;
104  }
105  }
106  work = r;
107  work *= 1/rnrm2;
108  V.push_back(work);
109  s.push_back(rnrm2);
110 
111  int jiter;
112  // FGMRES iteration
113  Log::info() << std::endl;
114  for (jiter = 0; jiter < maxiter; ++jiter) {
115  Log::info() << " FGMRES Starting Iteration " << jiter+1 << std::endl;
116 
117  precond.multiply(V[jiter], z);
118  Z.push_back(z);
119  A.multiply(z, w);
120 
121  double avnrm2 = sqrt(dot_product(w, w));
122 
123  // Modified Gram-Schmidt
124  for (int jj = 0; jj <= jiter; ++jj) {
125  H[jiter][jj] = dot_product(V[jj], w);
126  w.axpy(-H[jiter][jj], V[jj]);
127  }
128  H[jiter][jiter+1] = sqrt(dot_product(w, w));
129  double av2nrm2 = H[jiter][jiter+1];
130 
131  // Re-orthogonalize if necessary
132  if (avnrm2 + 0.001*av2nrm2 == avnrm2) {
133  for (int jj = 0; jj <= jiter; ++jj) {
134  double hr = dot_product(V[jj], w);
135  H[jiter][jj] += hr;
136  w.axpy(-hr, V[jj]);
137  }
138  H[jiter][jiter+1] = sqrt(dot_product(w, w));
139  }
140 
141  // Check breakdown
142  if (H[jiter][jiter+1] != 0) {
143  w *= 1/H[jiter][jiter+1];
144  }
145  V.push_back(w);
146 
147  if (jiter > 0) {
148  // Apply Givens Rotation
149  for (int jj = 0; jj < jiter; ++jj) {
150  double temp = cs[jj]*H[jiter][jj] + sn[jj]*H[jiter][jj+1];
151  H[jiter][jj+1] = -sn[jj]*H[jiter][jj] + cs[jj]*H[jiter][jj+1];
152  H[jiter][jj] = temp;
153  }
154  }
155 
156  // Compute Givens rotation matrix parameters
157  rotmat(H[jiter][jiter], H[jiter][jiter+1], cs[jiter], sn[jiter]);
158 
159  H[jiter][jiter] = cs[jiter]*H[jiter][jiter] + sn[jiter]*H[jiter][jiter+1];
160  H[jiter][jiter+1] = 0.0;
161 
162  double temp = cs[jiter]*s[jiter];
163  s.push_back(-sn[jiter]*s[jiter]);
164  s[jiter] = temp;
165 
166  // Calculate the solution
167  UpTriSolve(H, s, y, jiter); // H s = y
168  x = xinit;
169  for (int jj = 0; jj < jiter+1; ++jj) {
170  x.axpy(y[jj], Z[jj]);
171  }
172 
173  normReduction = std::abs(s[jiter+1])/bnrm2;
174  Log::info() << "FGMRES end of iteration " << jiter+1 << std::endl;
175  printNormReduction(jiter+1, std::abs(s[jiter+1]), normReduction);
176 
177  if (normReduction <= tolerance) {
178  Log::info() << "FGMRES: Achieved required reduction in residual norm." << std::endl;
179  jiter += 1;
180  break;
181  }
182  }
183 
184  // Calculate the solution
185  UpTriSolve(H, s, y, jiter); // H s = y
186  x = xinit;
187  for (int jj = 0; jj < jiter; ++jj) {
188  x.axpy(y[jj], Z[jj]);
189  }
190 
191  Log::info() << "FGMRES: end" << std::endl;
192 
193  return normReduction;
194 }
195 
196 } // namespace oops
197 
198 #endif // OOPS_ASSIMILATION_FGMRES_H_
Solves the upper triangular system sol = H \ rhs.
The namespace for the main oops code.
void UpTriSolve(const std::vector< std::vector< double > > &H, const std::vector< double > &rhs, std::vector< double > &sol, const int &dim)
Definition: UpTriSolve.h:22
void rotmat(const double &a, const double &b, double &c, double &s)
void printNormReduction(int iteration, const double &grad, const double &norm)
double FGMRES(VECTOR &x, const VECTOR &b, const AMATRIX &A, const PMATRIX &precond, const int maxiter, const double tolerance)
Definition: FGMRES.h:69
Compute the Givens rotation matrix parameters.