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