IODA Bundle
oops/assimilation/FullGMRES.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_FULLGMRES_H_
12 #define OOPS_ASSIMILATION_FULLGMRES_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 FullGMRES.h
28  * \brief FullGMRES solver for Ax=b.
29  *
30  * GMRES solver for Ax=b.(based on implementation following
31  * Saad-Schultz and C. T. Kelley, July 10, 1994)
32 
33  * A must be square, but need not be symmetric.
34  * A preconditioner must be supplied that, given a
35  * vector q, returns an approximate solution of Ap=q.
36  *
37  * On entry:
38  * - x = starting point, \f$ X_0 \f$.
39  * - b = right hand side.
40  * - A = \f$ A \f$.
41  * - precond = preconditioner \f$ F_k \approx (A)^{-1} \f$.
42  * - maxiter = maximum number of iterations
43  * - tol = error tolerance
44  *
45  * On exit, x will contain the solution.
46  * The return value is the achieved reduction in residual norm.
47  *
48  * Iteration will stop if the maximum iteration limit "maxiter" is reached
49  * or if the residual norm reduces by a factor of "tolerance".
50  *
51  * VECTOR must implement:
52  * - dot_product
53  * - operator(=)
54  * - operator(+=),
55  * - operator(-=)
56  * - operator(*=) [double * VECTOR],
57  * - axpy
58  *
59  * AMATRIX and PMATRIX must implement a method:
60  * - void multiply(const VECTOR&, VECTOR&) const
61  *
62  * which applies the matrix to the first argument, and returns the
63  * matrix-vector product in the second. (Note: the const is optional, but
64  * recommended.)
65  */
66 
67 template <typename VECTOR, typename AMATRIX, typename PMATRIX>
68 double FullGMRES(VECTOR & xx, const VECTOR & bb, const AMATRIX & A,
69  const PMATRIX & precond, const int maxiter,
70  const double tolerance, std::vector<VECTOR> & pqVEC,
71  std::vector<VECTOR> & xyVEC) {
72  std::vector<VECTOR> VV;
73  // reserve space to avoid extra copies
74  VV.reserve(maxiter+1);
75 
76  std::vector< std::vector<double> > H;
77  std::vector<double> cs(maxiter+1, 0);
78  std::vector<double> sn(maxiter+1, 0);
79  std::vector<double> ss;
80  std::vector<double> yy(maxiter+1, 0);
81 
82  VECTOR ww(xx);
83  VECTOR zz(xx);
84 
85  VECTOR rr(bb);
86  double xnrm2 = dot_product(xx, xx);
87  if (xnrm2 > 0.0) {
88  A.multiply(xx, ww);
89  rr -= ww; // r = b - Ax
90  }
91 
92  precond.multiply(rr, zz);
93 
94  double znrm2 = sqrt(dot_product(zz, zz));
95  double normReduction = 1.0;
96 
97  zz *= 1/znrm2;
98  VV.push_back(zz);
99  ss.push_back(znrm2);
100 
101  // Initialiaze (maxiter + 1) by maxiter matrix H
102  H.resize(maxiter);
103  for (int ii = 0; ii <= maxiter-1; ii++) {
104  H[ii].resize(maxiter + 1);
105  for (int jj = 0; jj <= maxiter; jj++) {
106  H[ii][jj] = 0;
107  }
108  }
109 
110  pqVEC.clear();
111  xyVEC.clear();
112 
113  int jiter;
114  // FullGMRES iteration
115  Log::info() << std::endl;
116  for (jiter = 0; jiter < maxiter; ++jiter) {
117  Log::info() << " FullGMRES Starting Iteration " << jiter+1 << std::endl;
118 
119  A.multiply(VV[jiter], zz);
120  precond.multiply(zz, ww);
121  if (jiter > 1) {
122  xyVEC.push_back(VV[jiter]);
123  pqVEC.push_back(zz);
124  }
125 
126  double avnrm2 = sqrt(dot_product(ww, ww));
127 
128  // Modified Gram-Schmidt
129  for (int jj = 0; jj <= jiter; ++jj) {
130  H[jiter][jj] = dot_product(VV[jj], ww);
131  ww.axpy(-H[jiter][jj], VV[jj]);
132  }
133 
134  H[jiter][jiter+1] = sqrt(dot_product(ww, ww));
135 
136  double av2nrm2 = H[jiter][jiter+1];
137 
138  // Re-orthogonalize if necessary
139  if (avnrm2 + 0.001*av2nrm2 == avnrm2) {
140  for (int jj = 0; jj <= jiter; ++jj) {
141  double hr = dot_product(VV[jj], ww);
142  H[jiter][jj] += hr;
143  ww.axpy(-hr, VV[jj]);
144  }
145  H[jiter][jiter+1] = sqrt(dot_product(ww, ww));
146  }
147 
148  // Check breakdown
149  if (H[jiter][jiter+1] != 0.0) {
150  ww *= 1/H[jiter][jiter+1];
151  }
152 
153  VV.push_back(ww);
154 
155  if (jiter > 0) {
156  // Apply Givens Rotation
157  for (int jj = 0; jj < jiter; ++jj) {
158  double temp = cs[jj]*H[jiter][jj] + sn[jj]*H[jiter][jj+1];
159  H[jiter][jj+1] = -sn[jj]*H[jiter][jj] + cs[jj]*H[jiter][jj+1];
160  H[jiter][jj] = temp;
161  }
162  }
163 
164  // Compute Givens rotation matrix parameters
165  rotmat(H[jiter][jiter], H[jiter][jiter+1], cs[jiter], sn[jiter]);
166 
167  H[jiter][jiter] = cs[jiter]*H[jiter][jiter] + sn[jiter]*H[jiter][jiter+1];
168  H[jiter][jiter+1] = 0.0;
169 
170  double temp = cs[jiter]*ss[jiter];
171  ss.push_back(-sn[jiter]*ss[jiter]);
172  ss[jiter] = temp;
173 
174  normReduction = std::abs(ss[jiter+1])/znrm2;
175  Log::info() << "FullGMRES end of iteration " << jiter+1 << std::endl;
176  printNormReduction(jiter+1, std::abs(ss[jiter+1]), normReduction);
177 
178  if (normReduction <= tolerance) {
179  Log::info() << "FullGMRES: Achieved required reduction in presidual norm." << std::endl;
180  jiter += 1;
181  break;
182  }
183  }
184 
185  // Calculate the solution
186  UpTriSolve(H, ss, yy, jiter);
187  for (int jj = 0; jj < jiter; ++jj) {
188  xx.axpy(yy[jj], VV[jj]);
189  }
190 
191  return normReduction;
192 }
193 
194 } // namespace oops
195 
196 #endif // OOPS_ASSIMILATION_FULLGMRES_H_
Solves the upper triangular system sol = H \ rhs.
The namespace for the main oops code.
double FullGMRES(VECTOR &xx, const VECTOR &bb, const AMATRIX &A, const PMATRIX &precond, const int maxiter, const double tolerance, std::vector< VECTOR > &pqVEC, std::vector< VECTOR > &xyVEC)
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)
Compute the Givens rotation matrix parameters.