OOPS
DRIPCGMinimizer.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_DRIPCGMINIMIZER_H_
12 #define OOPS_ASSIMILATION_DRIPCGMINIMIZER_H_
13 
14 #include <cmath>
15 #include <string>
16 #include <vector>
17 
25 #include "oops/util/dot_product.h"
26 #include "oops/util/formats.h"
27 #include "oops/util/Logger.h"
28 #include "oops/util/printRunStats.h"
29 
30 namespace oops {
31 
32 /// Derber-Rosati IPCG Minimizer
33 /*!
34  * \brief Derber-Rosati Inexact-Preconditioned Conjugate Gradients solver.
35  *
36  * This solver is based on the Golub-Ye Inexact-Preconditioned Conjugate
37  * Gradients solver for Ax=b (G.H. Golub and Q. Ye 1999/00,
38  * SIAM J. Sci. Comput. 21(4) 1305-1320), and on the Derber and Rosati double
39  * PCG algorithm (J. Derber and A. Rosati, 1989, J. Phys. Oceanog. 1333-1347).
40  * It solves \f$ Ax=b\f$ for the particular case \f$ A=B^{-1}+C\f$,
41  * without requiring the application of \f$ B^{-1}\f$.
42  *
43  * A must be square, symmetric, positive definite.
44  *
45  * A preconditioner must be supplied that, given a vector q, returns an
46  * approximation to \f$ (AB)^{-1} q\f$. The preconditioner can be variable.
47  * Note that the traditional \f$ B\f$-preconditioning corresponds to
48  * precond=\f$I\f$.
49  *
50  * On entry:
51  * - x = starting point, \f$ X_0 \f$.
52  * - xh = \f$ B^{-1} x_0\f$.
53  * - b = right hand side.
54  * - B = \f$ B \f$.
55  * - C = \f$ C \f$.
56  * - precond = preconditioner \f$ F_k \approx (AB)^{-1} \f$.
57  *
58  * On exit, x will contain the solution \f$ x \f$ and xh will contain
59  * \f$ B^{-1} x\f$.
60  * The return value is the achieved reduction in residual norm.
61  *
62  * Iteration will stop if the maximum iteration limit "maxiter" is reached
63  * or if the residual norm reduces by a factor of "tolerance".
64  *
65  * Each of BMATRIX, CMATRIX and PMATRIX must implement a method:
66  * - void multiply(const VECTOR&, VECTOR&) const
67  *
68  * which applies the matrix to the first argument, and returns the
69  * matrix-vector product in the second. (Note: the const is optonal, but
70  * recommended.)
71  */
72 
73 // -----------------------------------------------------------------------------
74 
75 template<typename MODEL, typename OBS> class DRIPCGMinimizer : public DRMinimizer<MODEL, OBS> {
80 
81  public:
82  const std::string classname() const override {return "DRIPCGMinimizer";}
83  DRIPCGMinimizer(const eckit::Configuration &, const CostFct_ &);
85 
86  private:
87  double solve(CtrlInc_ &, CtrlInc_ &, CtrlInc_ &, const Bmat_ &, const HtRinvH_ &,
88  const double, const double, const int, const double) override;
89 
91 };
92 
93 // =============================================================================
94 
95 template<typename MODEL, typename OBS>
96 DRIPCGMinimizer<MODEL, OBS>::DRIPCGMinimizer(const eckit::Configuration & conf, const CostFct_ & J)
97  : DRMinimizer<MODEL, OBS>(J), lmp_(conf)
98 {}
99 
100 // -----------------------------------------------------------------------------
101 
102 template<typename MODEL, typename OBS>
104  const Bmat_ & B, const HtRinvH_ & HtRinvH,
105  const double costJ0Jb, const double costJ0JoJc,
106  const int maxiter, const double tolerance) {
107  util::printRunStats("DRIPCG start");
108  CtrlInc_ ap(xh);
109  CtrlInc_ pp(xh);
110  CtrlInc_ ph(xh);
111  CtrlInc_ ss(xh);
112  CtrlInc_ sh(xh);
113  CtrlInc_ dr(xh);
114  CtrlInc_ r0(xh);
115 
116  std::vector<CtrlInc_> vvecs; // for re-orthogonalization
117  std::vector<CtrlInc_> zvecs; // for re-orthogonalization
118  std::vector<double> scals; // for re-orthogonalization
119  // reserve space in vectors to avoid extra copies
120  vvecs.reserve(maxiter+1);
121  zvecs.reserve(maxiter+1);
122  scals.reserve(maxiter+1);
123 
124  const double costJ0 = costJ0Jb + costJ0JoJc;
125 
126  r0 = rr;
127 
128  lmp_.multiply(rr, sh);
129  B.multiply(sh, ss);
130 
131  double rrnorm0 = sqrt(dot_product(rr, rr));
132  double dotSr0 = dot_product(rr, ss);
133  double normReduction = 1.0;
134  double rdots = dotSr0;
135  double rdots_old = dotSr0;
136 
137  vvecs.push_back(rr);
138  zvecs.push_back(ss);
139  scals.push_back(1.0/dotSr0);
140 
141  Log::info() << std::endl;
142  for (int jiter = 0; jiter < maxiter; ++jiter) {
143  Log::info() << " DRIPCG Starting Iteration " << jiter+1 << std::endl;
144  util::printRunStats("DRIPCG iteration " + std::to_string(jiter+1));
145 
146  if (jiter == 0) {
147  pp = ss;
148  ph = sh;
149  } else {
150  dr -= rr; // dr=oldr-r
151  double beta = -dot_product(ss, dr)/rdots_old;
152 
153  pp *= beta;
154  pp += ss; // p = s + beta*p
155 
156  ph *= beta;
157  ph += sh; // ph = sh + beta*ph
158  }
159 
160  HtRinvH.multiply(pp, ap);
161  ap += ph;
162 
163  dr = rr;
164 
165  double rho = dot_product(pp, ap);
166  double alpha = rdots/rho;
167  Log::info() << "DRIPCGMinimizer rho = " << rho << ", alpha = " << alpha << std::endl;
168 
169  xx.axpy(alpha, pp); // xx = xx + alpha*pp
170  xh.axpy(alpha, ph); // xh = xh + alpha*ph
171  rr.axpy(-alpha, ap); // rr = rr - alpha*ap
172 
173  // Compute the quadratic cost function
174  double costJ = costJ0 - 0.5 * dot_product(xx, r0);
175  double costJb = costJ0Jb + 0.5 * dot_product(xx, xh);
176  double costJoJc = costJ - costJb;
177 
178  // Re-orthogonalization
179  for (int jj = 0; jj < jiter; ++jj) {
180  double proj = scals[jj] * dot_product(rr, zvecs[jj]);
181  rr.axpy(-proj, vvecs[jj]);
182  }
183 
184  lmp_.multiply(rr, sh);
185  B.multiply(sh, ss);
186 
187  rdots_old = rdots;
188  rdots = dot_product(rr, ss);
189 // There is an issue where in some cases ss goes to zero before rr. The convergence
190 // check below only checks for rr. Leaving the commented lines here for further invertigation.
191 // double sdots = dot_product(ss, ss);
192 // double rdotr = dot_product(rr, rr);
193 // Log::info() << "DRIPCGMinimizer rdots = " << rdots
194 // << ", sdots = " << sdots << ", rdotr = " << rdotr << std::endl;
195  double rrnorm = sqrt(dot_product(rr, rr));
196  normReduction = rrnorm/rrnorm0;
197 
198  Log::info() << "DRIPCG end of iteration " << jiter+1 << std::endl;
199  printNormReduction(jiter+1, rrnorm, normReduction);
200  printQuadraticCostFunction(jiter+1, costJ, costJb, costJoJc);
201 
202  // Save the pairs for preconditioning
203  lmp_.push(pp, ph, ap, rho);
204 
205  if (normReduction < tolerance) {
206  Log::info() << "DRIPCG: Achieved required reduction in residual norm." << std::endl;
207  break;
208  }
209 
210  vvecs.push_back(rr);
211  zvecs.push_back(ss);
212  scals.push_back(1.0/rdots);
213  }
214 
215 // Generate the (second-level) Limited Memory Preconditioner
216  lmp_.update(B);
217 
218  util::printRunStats("DRIPCG end");
219  return normReduction;
220 }
221 
222 // -----------------------------------------------------------------------------
223 
224 } // namespace oops
225 
226 #endif // OOPS_ASSIMILATION_DRIPCGMINIMIZER_H_
The matrix.
Definition: BMatrix.h:27
void multiply(const CtrlInc_ &x, CtrlInc_ &bx) const
Definition: BMatrix.h:33
void axpy(const double, const ControlIncrement &)
Cost Function.
Definition: CostFunction.h:53
Derber-Rosati IPCG Minimizer.
DRIPCGMinimizer(const eckit::Configuration &, const CostFct_ &)
BMatrix< MODEL, OBS > Bmat_
HtRinvHMatrix< MODEL, OBS > HtRinvH_
CostFunction< MODEL, OBS > CostFct_
const std::string classname() const override
ControlIncrement< MODEL, OBS > CtrlInc_
double solve(CtrlInc_ &, CtrlInc_ &, CtrlInc_ &, const Bmat_ &, const HtRinvH_ &, const double, const double, const int, const double) override
QNewtonLMP< CtrlInc_, Bmat_ > lmp_
DR (Derber and Rosati) Minimizers.
Definition: DRMinimizer.h:46
void multiply(const CtrlInc_ &dx, CtrlInc_ &dz) const
Definition: HtRinvHMatrix.h:63
The namespace for the main oops code.
void printQuadraticCostFunction(int iteration, const double &costJ, const double &costJb, const double &costJoJc)
void printNormReduction(int iteration, const double &grad, const double &norm)