Loading [MathJax]/extensions/tex2jax.js
OOPS
All Classes Namespaces Files Functions Variables Typedefs Macros Pages
PCG.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_PCG_H_
12 #define OOPS_ASSIMILATION_PCG_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 PCG.h
25  * \brief Preconditioned Conjugate Gradients solver.
26  *
27  * This solver is based on the standard Preconditioned Conjugate
28  * Gradients solver for Ax=b
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.
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$ P \approx (A)^{-1} \f$.
39  *
40  * On exit, x will contain the solution \f$ x \f$
41 
42  * The return value is the achieved reduction in preconditioned residual norm.
43  *
44  * Iteration will stop if the maximum iteration limit "maxiter" is reached
45  * or if the preconditioned 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  * Each of 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 PCG(VECTOR & x, const VECTOR & b,
65  const AMATRIX & A, const PMATRIX & precond,
66  const int maxiter, const double tolerance ) {
67  VECTOR ap(x);
68  VECTOR p(x);
69  VECTOR r(x);
70  VECTOR s(x);
71  VECTOR v(x);
72  VECTOR z(x);
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, s);
92  double normReduction = 1.0;
93  double rdots = dotRr0;
94  double rdots_old = 0.0;
95 
96  v = r;
97  v *= 1/sqrt(dotRr0);
98  z = s;
99  z *= 1/sqrt(dotRr0);
100 
101  vVEC.push_back(v);
102  zVEC.push_back(z);
103 
104  Log::info() << std::endl;
105  for (int jiter = 0; jiter < maxiter; ++jiter) {
106  Log::info() << " PCG Starting Iteration " << jiter+1 << std::endl;
107 
108  if (jiter == 0) {
109  p = s;
110  } else {
111  double beta = dot_product(s, r)/rdots_old;
112  Log::info() << "PCG beta = " << beta << std::endl;
113 
114  p *= beta;
115  p += s; // p = s + beta*p
116  }
117 
118  A.multiply(p, ap); // ap = A*p
119 
120  double alpha = rdots/dot_product(p, ap);
121 
122  x.axpy(alpha, p); // x = x + alpha*p
123  r.axpy(-alpha, ap); // r = r - alpha*ap
124 
125  // Re-orthogonalization
126  for (int iiter = 0; iiter < jiter; ++iiter) {
127  double proj = dot_product(r, zVEC[iiter]);
128  r.axpy(-proj, vVEC[iiter]);
129  }
130 
131  precond.multiply(r, s);
132 
133  rdots_old = rdots;
134  rdots = dot_product(r, s);
135 
136  v = r;
137  v *= 1/sqrt(rdots);
138  z = s;
139  z *= 1/sqrt(rdots);
140  vVEC.push_back(v);
141  zVEC.push_back(z);
142 
143  normReduction = sqrt(rdots/dotRr0);
144 
145  Log::info() << "PCG end of iteration " << jiter+1 << std::endl;
146  printNormReduction(jiter+1, sqrt(rdots), normReduction);
147 
148  if (normReduction < tolerance) {
149  Log::info() << "PCG: Achieved required reduction in residual norm." << std::endl;
150  break;
151  }
152  }
153 
154  Log::info() << "PCG: end" << std::endl;
155 
156  return normReduction;
157 }
158 
159 } // namespace oops
160 
161 #endif // OOPS_ASSIMILATION_PCG_H_
The namespace for the main oops code.
double PCG(VECTOR &x, const VECTOR &b, const AMATRIX &A, const PMATRIX &precond, const int maxiter, const double tolerance)
Definition: PCG.h:64
void printNormReduction(int iteration, const double &grad, const double &norm)