OOPS
IPCG.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_IPCG_H_
12 #define OOPS_ASSIMILATION_IPCG_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 IPCG.h
25  * \brief Inexact-Preconditioned Conjugate Gradients solver.
26  *
27  * Golub-Ye Inexact-Preconditioned Conjugate Gradients solver for Ax=b.
28  * (G.H. Golub and Q. Ye 1999/00, SIAM J. Sci. Comput. 21(4) 1305-1320.)
29  *
30  * A must be square, symmetric, positive definite.
31  * A preconditioner must be supplied that, given a vector q, returns an
32  * approximate solution of Ap=q. The preconditioner can be variable.
33  *
34  * On entry:
35  * - x = starting point, \f$ x_0 \f$.
36  * - b = right hand side.
37  * - A = \f$ A \f$.
38  * - precond = preconditioner \f$ F_k \approx (A)^{-1} \f$.
39  *
40  * On exit, x will contain the solution.
41  * The return value is the achieved reduction in residual norm.
42  *
43  * Iteration will stop if the maximum iteration limit "maxiter" is reached
44  * or if the residual norm reduces by a factor of "tolerance".
45  *
46  * VECTOR must implement:
47  * - dot_product
48  * - operator(=)
49  * - operator(+=),
50  * - operator(-=)
51  * - operator(*=) [double * VECTOR],
52  * - axpy
53  *
54  * AMATRIX and PMATRIX must implement a method:
55  * - void multiply(const VECTOR&, VECTOR&) const
56  *
57  * which applies the matrix to the first argument, and returns the
58  * matrix-vector product in the second. (Note: the const is optional, but
59  * recommended.)
60  */
61 
62 template <typename VECTOR, typename AMATRIX, typename PMATRIX>
63 double IPCG(VECTOR & x, const VECTOR & b,
64  const AMATRIX & A, const PMATRIX & precond,
65  const int maxiter, const double tolerance ) {
66  VECTOR ap(x);
67  VECTOR p(x);
68  VECTOR r(x);
69  VECTOR s(x);
70  VECTOR oldr(x);
71  VECTOR w(x);
72  VECTOR v(x); // required for re-orthogonalization
73  VECTOR z(x); // required for re-orthogonalization
74 
75  std::vector<VECTOR> vVEC; // required for re-orthogonalization
76  std::vector<VECTOR> zVEC; // required for re-orthogonalization
77  // reserve space to avoid extra copies
78  vVEC.reserve(maxiter+1);
79  zVEC.reserve(maxiter+1);
80 
81  // Initial residual r = b - Ax
82  r = b;
83  double xnrm2 = dot_product(x, x);
84  if (xnrm2 > 0.0) {
85  A.multiply(x, s);
86  r -= s;
87  }
88 
89  // s = precond r
90  precond.multiply(r, s);
91 
92  double rnorm0 = sqrt(dot_product(r, r));
93  double dotSr0 = dot_product(r, s);
94  double normReduction = 1.0;
95  double rdots_old = dotSr0;
96  double rdots = dotSr0;
97 
98  v = r;
99  v *= 1/sqrt(dotSr0);
100  z = s;
101  z *= 1/sqrt(dotSr0);
102 
103  vVEC.push_back(v);
104  zVEC.push_back(z);
105 
106  Log::info() << std::endl;
107  for (int jiter = 0; jiter < maxiter; ++jiter) {
108  Log::info() << " IPCG Starting Iteration " << jiter+1 << std::endl;
109 
110  if (jiter == 0) {
111  p = s;
112  } else {
113  w = r;
114  w -= oldr; // w=r-oldr
115  double beta = dot_product(s, w)/rdots_old;
116  p *= beta;
117  p += s; // p = s + beta*p
118  }
119 
120  A.multiply(p, ap); // ap = Ap
121 
122  oldr = r;
123 
124  double alpha = rdots/dot_product(p, ap);
125 
126  x.axpy(alpha, p); // x = x + alpha*p;
127  r.axpy(-alpha, ap); // r = r - alpha*ap;
128 
129  // Re-orthogonalization
130  for (int iiter = 0; iiter < jiter; ++iiter) {
131  double proj = dot_product(r, zVEC[iiter]);
132  r.axpy(-proj, vVEC[iiter]);
133  }
134 
135  precond.multiply(r, s); // returns s as approximate solve of As=r
136 
137  rdots_old = rdots;
138  rdots = dot_product(r, s);
139 
140  v = r;
141  v *= 1/sqrt(rdots);
142  z = s;
143  z *= 1/sqrt(rdots);
144  vVEC.push_back(v);
145  zVEC.push_back(z);
146 
147  double rnorm = sqrt(dot_product(r, r));
148  normReduction = rnorm/rnorm0;
149  Log::info() << "IPCG end of iteration " << jiter+1 << std::endl;
150  printNormReduction(jiter+1, rnorm, normReduction);
151 
152  if (normReduction < tolerance) {
153  Log::info() << "IPCG: Achieved required reduction in residual norm." << std::endl;
154  break;
155  }
156  }
157 
158  Log::info() << "IPCG: end" << std::endl;
159 
160  return normReduction;
161 }
162 
163 } // namespace oops
164 
165 #endif // OOPS_ASSIMILATION_IPCG_H_
The namespace for the main oops code.
double IPCG(VECTOR &x, const VECTOR &b, const AMATRIX &A, const PMATRIX &precond, const int maxiter, const double tolerance)
Definition: IPCG.h:63
void printNormReduction(int iteration, const double &grad, const double &norm)