OOPS
PLanczos.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_PLANCZOS_H_
12 #define OOPS_ASSIMILATION_PLANCZOS_H_
13 
14 #include <cmath>
15 #include <vector>
16 
19 #include "oops/util/dot_product.h"
20 #include "oops/util/formats.h"
21 #include "oops/util/Logger.h"
22 
23 namespace oops {
24 
25 /*! \file PLanczos.h
26  * \brief Preconditioned Lanczos solver.
27  *
28  * This algorihtm is based on the standard Preconditioned Lanczos
29  * method for solving the linear system Ax = b
30  *
31  * A must be square and symmetric.
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$ P_k \f$.
38  *
39  * On exit, x will contain the solution \f$ x \f$
40  *
41  * The return value is the achieved reduction in preconditioned 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  * Each of 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 PLanczos(VECTOR & xx, const VECTOR & bb,
64  const AMATRIX & A, const PMATRIX & precond,
65  const int maxiter, const double tolerance) {
66  VECTOR zz(xx);
67  VECTOR ww(xx);
68 
69  std::vector<VECTOR> vVEC;
70  std::vector<VECTOR> zVEC;
71  // reserve space to avoid extra copies
72  vVEC.reserve(maxiter+1);
73  zVEC.reserve(maxiter+1);
74 
75  Log::debug() << "PLanczos xx = " << xx << std::endl;
76  Log::debug() << "PLanczos bb = " << bb << std::endl;
77 
78  std::vector<double> alphas;
79  std::vector<double> betas;
80  std::vector<double> dd;
81  std::vector<double> yy;
82 
83  // Initial residual r = b - Ax
84  VECTOR rr(bb);
85  double xnrm2 = dot_product(xx, xx);
86  if (xnrm2 > 0.0) {
87  A.multiply(xx, zz);
88  rr -= zz;
89  }
90  Log::debug() << "PLanczos rr = " << rr << std::endl;
91 
92  // z = precond r
93  precond.multiply(rr, zz);
94  Log::debug() << "PLanczos zz = " << zz << std::endl;
95 
96  double normReduction = 1.0;
97  double beta0 = sqrt(dot_product(rr, zz));
98  double beta = 0.0;
99  Log::debug() << "PLanczos beta0 = " << beta0 << std::endl;
100 
101  VECTOR vv(rr);
102  vv *= 1/beta0;
103  zz *= 1/beta0;
104 
105  zVEC.push_back(zz); // zVEC[0] = z_1 ---> required for re-orthogonalization
106  vVEC.push_back(vv); // vVEC[0] = v_1 ---> required for re-orthogonalization
107 
108  Log::debug() << "PLanczos vv = " << vv << std::endl;
109  Log::debug() << "PLanczos zz = " << zz << std::endl;
110 
111  int jiter = 0;
112  Log::info() << std::endl;
113  while (jiter < maxiter) {
114  Log::info() << " PLanczos Starting Iteration " << jiter+1 << std::endl;
115 
116  // w = A z - beta * vold
117  A.multiply(zz, ww); // w = A z
118  Log::debug() << "PLanczos ww = " << ww << std::endl;
119  if (jiter > 0) ww.axpy(-beta, vVEC[jiter-1]);
120 
121  double alpha = dot_product(zz, ww);
122  Log::debug() << "PLanczos alpha = " << alpha << std::endl;
123 
124  ww.axpy(-alpha, vv); // w = w - alpha * v
125  Log::debug() << "PLanczos ww = " << ww << std::endl;
126 
127  // Re-orthogonalization
128  for (int iiter = 0; iiter < jiter; ++iiter) {
129  double proj = dot_product(ww, zVEC[iiter]);
130  ww.axpy(-proj, vVEC[iiter]);
131  }
132  Log::debug() << "PLanczos ww = " << ww << std::endl;
133 
134  precond.multiply(ww, zz); // z = precond w
135  Log::debug() << "PLanczos zz = " << zz << std::endl;
136 
137  beta = sqrt(dot_product(zz, ww));
138  Log::debug() << "PLanczos beta = " << beta << std::endl;
139 
140  vv = ww;
141  vv *= 1/beta;
142  zz *= 1/beta;
143  Log::debug() << "PLanczos vv = " << vv << std::endl;
144  Log::debug() << "PLanczos zz = " << zz << std::endl;
145 
146  zVEC.push_back(zz); // zVEC[jiter+1] = z_jiter
147  vVEC.push_back(vv); // vVEC[jiter+1] = v_jiter
148 
149  alphas.push_back(alpha);
150 
151  if (jiter == 0) {
152  yy.push_back(beta0/alpha);
153  dd.push_back(beta0);
154  } else {
155  // Solve the tridiagonal system T_jiter y_jiter = beta0 * e_1
156  dd.push_back(beta0*dot_product(zVEC[0], vv));
157  TriDiagSolve(alphas, betas, dd, yy);
158  }
159 
160  // Gradient norm in precond metric --> sqrt(r'z) --> beta * y(jiter)
161  double rznorm = beta*std::abs(yy[jiter]);
162  Log::debug() << "PLanczos rznorm = " << rznorm << std::endl;
163 
164  normReduction = rznorm/beta0;
165 
166  betas.push_back(beta);
167 
168  Log::info() << "PLanczos end of iteration " << jiter+1 << std::endl;
169  printNormReduction(jiter+1, rznorm, normReduction);
170 
171  ++jiter;
172 
173  if (normReduction < tolerance) {
174  Log::info() << "PLanczos: Achieved required reduction in residual norm." << std::endl;
175  break;
176  }
177  }
178  Log::debug() << "PLanczos jiter = " << jiter << std::endl;
179 
180  // Calculate the solution (xh = Binv x)
181  for (int iiter = 0; iiter < jiter; ++iiter) {
182  xx.axpy(yy[iiter], zVEC[iiter]);
183  }
184 
185  Log::info() << "PLanczos: end" << std::endl;
186 
187  return normReduction;
188 }
189 
190 } // namespace oops
191 
192 #endif // OOPS_ASSIMILATION_PLANCZOS_H_
The namespace for the main oops code.
void TriDiagSolve(const std::vector< double > &diag, const std::vector< double > &sub, const std::vector< double > &rhs, std::vector< double > &sol)
double PLanczos(VECTOR &xx, const VECTOR &bb, const AMATRIX &A, const PMATRIX &precond, const int maxiter, const double tolerance)
Definition: PLanczos.h:63
void printNormReduction(int iteration, const double &grad, const double &norm)