OOPS
GMRESR.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_GMRESR_H_
12 #define OOPS_ASSIMILATION_GMRESR_H_
13 
14 #include <cmath>
15 #include <vector>
16 
18 #include "oops/util/dot_product.h"
19 #include "oops/util/formats.h"
20 #include "oops/util/Logger.h"
21 
22 namespace oops {
23 
24 /*! \file GMRESR.h
25  * \brief GMRESR solver for Ax=b.
26  *
27  * GMRESR solver for Ax=b. (H.A. Van der Vorst and C. Vuik, 1994,
28  * Numerical Linear Algebra with Applications, 1(4), 369-386.)
29  * A must be square, but need not be symmetric. (For a symmetric matrix,
30  * and constant preconditioner, GMRESR is simply PCG with full
31  * orthogonalisation.) A preconditioner must be supplied that, given a
32  * vector q, returns an approximate solution of Ap=q. The preconditioner
33  * can be variable.
34  *
35  * On entry:
36  * - x = starting point, \f$ X_0 \f$.
37  * - b = right hand side.
38  * - A = \f$ A \f$.
39  * - precond = preconditioner \f$ F_k \approx (A)^{-1} \f$.
40  *
41  * On exit, x will contain the solution.
42  * The return value is the achieved reduction in residual norm.
43  *
44  * Iteration will stop if the maximum iteration limit "maxiter" is reached
45  * or if the residual norm reduces by a factor of "tolerance".
46  *
47  * VECTOR must implement:
48  * - dot_product
49  * - operator(=)
50  * - operator(+=),
51  * - operator(-=)
52  * - operator(*=) [double * VECTOR],
53  * - axpy
54  *
55  * AMATRIX and PMATRIX must implement a method:
56  * - void multiply(const VECTOR&, VECTOR&) const
57  *
58  * which applies the matrix to the first argument, and returns the
59  * matrix-vector product in the second. (Note: the const is optional, but
60  * recommended.)
61  */
62 
63 template <typename VECTOR, typename AMATRIX, typename PMATRIX>
64 double GMRESR(VECTOR & xx, const VECTOR & bb,
65  const AMATRIX & A, const PMATRIX & precond,
66  const int maxiter, const double tolerance) {
67  const double stagthresh = 1.0e-3;
68  const double smallres = 1.0e-6;
69 
70  std::vector<VECTOR> cvecs;
71  std::vector<VECTOR> zvecs;
72  VECTOR cc(xx);
73  VECTOR zz(xx);
74 
75  VECTOR rr(bb);
76  A.multiply(xx, zz); // zz=Axx
77  rr -= zz;
78 
79  double normReduction = 1.0;
80  double rrnorm = sqrt(dot_product(rr, rr));
81  double cdotr = rrnorm;
82  const double rrnorm0 = rrnorm;
83 
84  if (rrnorm > smallres) {
85  Log::info() << std::endl;
86  for (int jiter = 0; jiter < maxiter; ++jiter) {
87  Log::info() << " GMRESR Starting Iteration " << jiter+1 << std::endl;
88 
89  // test for stagnation
90  if (std::abs(cdotr) < stagthresh*rrnorm) {
91  Log::info() << "GMRESR stagnated. Doing an LSQR step." << std::endl;
92  A.multiply(rr, zz); // should be A^T, but we assume A is symmetric.
93  } else {
94  precond.multiply(rr, zz); // returns zz as approximate solution of Azz=rr
95  }
96 
97  A.multiply(zz, cc); // cc=Azz
98 
99  for (int jj = 0; jj < jiter; ++jj) {
100  double alpha = -dot_product(cvecs[jj], cc);
101  cc.axpy(alpha, cvecs[jj]); // cc = cc - <c[jj],cc>*c[jj];
102  zz.axpy(alpha, zvecs[jj]); // zz = zz - <c[jj],cc>*z[jj];
103  }
104 
105  double ccnorm = sqrt(dot_product(cc, cc));
106  cvecs.push_back(cc); // c[jiter]=cc
107  cvecs[jiter] *= 1.0/ccnorm;
108  zvecs.push_back(zz); // z[jiter]=zz
109  zvecs[jiter] *= 1.0/ccnorm;
110 
111  cdotr = dot_product(cvecs[jiter], rr);
112  rr.axpy(-cdotr, cvecs[jiter]);
113  xx.axpy(cdotr, zvecs[jiter]);
114 
115  rrnorm = sqrt(dot_product(rr, rr));
116  normReduction = rrnorm/rrnorm0;
117  Log::info() << "GMRESR end of iteration " << jiter+1 << std::endl;
118  printNormReduction(jiter+1, rrnorm, normReduction);
119 
120  if (normReduction < tolerance) {
121  Log::info() << "GMRESR: Achieved required reduction in residual norm." << std::endl;
122  break;
123  }
124  }
125  }
126 
127  return normReduction;
128 }
129 
130 } // namespace oops
131 
132 #endif // OOPS_ASSIMILATION_GMRESR_H_
The namespace for the main oops code.
double GMRESR(VECTOR &xx, const VECTOR &bb, const AMATRIX &A, const PMATRIX &precond, const int maxiter, const double tolerance)
Definition: GMRESR.h:64
void printNormReduction(int iteration, const double &grad, const double &norm)