OOPS
DRPCGMinimizer.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_DRPCGMINIMIZER_H_
12 #define OOPS_ASSIMILATION_DRPCGMINIMIZER_H_
13 
14 #include <cmath>
15 #include <string>
16 #include <vector>
17 
18 #include "eckit/config/Configuration.h"
26 #include "oops/util/dot_product.h"
27 #include "oops/util/formats.h"
28 #include "oops/util/Logger.h"
29 
30 namespace oops {
31 
32 /// DRPCG Minimizer
33 /*!
34  * \brief Derber-Rosati Preconditioned Conjugate Gradients solver.
35  *
36  * This solver is based on the standard Preconditioned Conjugate
37  * Gradients solver for Ax=b (G. H. Golub and C. F. Van Loan, Matrix
38  * Computations), and on the Derber and Rosati double
39  * PCG algorithm (J. Derber and A. Rosati, 1989, J. Phys. Oceanog. 1333-1347).
40  * For details see S. Gurol, PhD Manuscript, 2013.
41  * It solves \f$ Ax=b\f$ for the particular case \f$ A=B^{-1}+C\f$,
42  * without requiring the application of \f$ B^{-1}\f$. This algorithm
43  * is similar to DRIPCG except it includes standard PCG instead IPCG
44  * and stopping criteria is based on the preconditioner norm.
45  *
46  * A must be square, symmetric, positive definite.
47  *
48  * A preconditioner must be supplied that, given a vector q, returns an
49  * approximation to \f$ (AB)^{-1} q\f$. Possible preconditioning
50  * is detailed in S. Gurol, PhD Manuscript, 2013.
51  * Note that the traditional \f$ B\f$-preconditioning corresponds to
52  * precond=\f$I\f$.
53  *
54  * On entry:
55  * - dx = starting point, \f$ dx_{0} \f$.
56  * - dxh = \f$ B^{-1} dx_{0} \f$.
57  * - rr = \f$ (sum B^-1 dx^{b}_{i} + ) H^T R^{-1} d \f$
58  * - B = \f$ B \f$.
59  * - C = \f$ C \f$.
60  * - precond = preconditioner \f$ F_k \approx (AB)^{-1} \f$.
61  *
62  * On exit, dx will contain the solution \f$ dx \f$ and dxh will contain
63  * \f$ B^{-1} dx\f$.
64  * The return value is the achieved reduction in preconditioned residual norm.
65  *
66  * Iteration will stop if the maximum iteration limit "maxiter" is reached
67  * or if the residual norm reduces by a factor of "tolerance".
68  *
69  * Each matrix must implement a method:
70  * - void multiply(const VECTOR&, VECTOR&) const
71  *
72  * which applies the matrix to the first argument, and returns the
73  * matrix-vector product in the second. (Note: the const is optional, but
74  * recommended.)
75  */
76 
77 // -----------------------------------------------------------------------------
78 
79 template<typename MODEL, typename OBS> class DRPCGMinimizer : public DRMinimizer<MODEL, OBS> {
84 
85  public:
86  const std::string classname() const override {return "DRPCGMinimizer";}
87  DRPCGMinimizer(const eckit::Configuration &, const CostFct_ &);
89 
90  private:
91  double solve(CtrlInc_ &, CtrlInc_ &, CtrlInc_ &, const Bmat_ &, const HtRinvH_ &,
92  const double, const double, const int, const double) override;
93 
95 };
96 
97 // =============================================================================
98 
99 template<typename MODEL, typename OBS>
100 DRPCGMinimizer<MODEL, OBS>::DRPCGMinimizer(const eckit::Configuration & conf, const CostFct_ & J)
101  : DRMinimizer<MODEL, OBS>(J), lmp_(conf)
102 {}
103 
104 // -----------------------------------------------------------------------------
105 
106 template<typename MODEL, typename OBS>
108  const Bmat_ & B, const HtRinvH_ & HtRinvH,
109  const double costJ0Jb, const double costJ0JoJc,
110  const int maxiter, const double tolerance) {
111  // dx increment
112  // dxh B^{-1} dx
113  // rr (sum B^{-1} dx_i^{b} +) G^T H^{-1} d
114 
115  CtrlInc_ qq(dxh);
116  CtrlInc_ pp(dxh);
117  CtrlInc_ hh(dxh);
118  CtrlInc_ zz(dxh);
119  CtrlInc_ pr(dxh);
120  CtrlInc_ r0(dxh);
121  CtrlInc_ ww(dxh);
122 
123  // vectors for re-orthogonalization
124  std::vector<CtrlInc_> vvecs;
125  std::vector<CtrlInc_> zvecs;
126  std::vector<double> scals;
127  // reserve space in vectors to avoid extra copies
128  vvecs.reserve(maxiter+1);
129  zvecs.reserve(maxiter+1);
130  scals.reserve(maxiter+1);
131 
132  // J0
133  const double costJ0 = costJ0Jb + costJ0JoJc;
134 
135  // r_{0}
136  r0 = rr;
137 
138  // r_{0}^T r_{0}
139  Log::info() << "normr0 " << dot_product(rr, rr) << std::endl;
140 
141  // z_{0} = B LMP r_{0}
142  lmp_.multiply(rr, pr);
143  B.multiply(pr, zz);
144  // p_{0} = z_{0}
145  pp = zz;
146  // h_{0} = LMP r_{0}
147  hh = pr;
148 
149  // r_{0}^T z_{0}
150  double dotRr0 = dot_product(rr, zz);
151  double normReduction = 1.0;
152  double rdots = dotRr0;
153  double rdots_old = 0.0;
154 
155  vvecs.push_back(rr);
156  zvecs.push_back(zz);
157  scals.push_back(1.0/dotRr0);
158 
159  Log::info() << std::endl;
160  for (int jiter = 0; jiter < maxiter; ++jiter) {
161  Log::info() << " DRPCG Starting Iteration " << jiter+1 << std::endl;
162 
163  if (jiter > 0) {
164  // beta_{i} = r_{i+1}^T z_{i+1} / r_{i}^T z_{i}
165  double beta = rdots/rdots_old;
166 
167  // p_{i+1} = z_{i+1} + beta*p_{i}
168  pp *= beta;
169  pp += zz;
170 
171  // h_{i+1} = LMP r_{i+1} + beta*h_{i}
172  hh *= beta;
173  hh += pr;
174  }
175 
176  // q_{i} = h_{i} + H^T R^{-1} H p_{i}
177  HtRinvH.multiply(pp, qq);
178  qq += hh;
179 
180  // alpha_{i} = r_{i}^T z_{i} / q_{i}^T p_{i}
181  double rho = dot_product(pp, qq);
182  double alpha = rdots/rho;
183 
184  // dx_{i+1} = dx_{i} + alpha * p_{i}
185  dx.axpy(alpha, pp);
186  // dxh_{i+1} = dxh_{i} + alpha * h_{i} ! for diagnosing Jb
187  dxh.axpy(alpha, hh);
188  // r_{i+1} = r_{i} - alpha * q_{i}
189  rr.axpy(-alpha, qq);
190 
191  // Compute the quadratic cost function
192  // J[dx_{i}] = J[0] - 0.5 dx_{i}^T r_{0}
193  double costJ = costJ0 - 0.5 * dot_product(dx, r0);
194  // Jb[dx_{i}] = 0.5 dx_{i}^T f_{i}
195  double costJb = costJ0Jb + 0.5 * dot_product(dx, dxh);
196  // Jo[dx_{i}] + Jc[dx_{i}] = J[dx_{i}] - Jb[dx_{i}]
197  double costJoJc = costJ - costJb;
198 
199  // Re-orthogonalization
200  for (int jj = 0; jj < jiter; ++jj) {
201  double proj = scals[jj] * dot_product(rr, zvecs[jj]);
202  rr.axpy(-proj, vvecs[jj]);
203  }
204 
205  // z_{i+1} = B LMP r_{i+1}
206  lmp_.multiply(rr, pr);
207  B.multiply(pr, zz);
208 
209  // r_{i}^T z_{i}
210  rdots_old = rdots;
211  // r_{i+1}^T z_{i+1}
212  rdots = dot_product(rr, zz);
213 
214  Log::info() << "rdots " << rdots << " iteration " << jiter << std::endl;
215 
216  // r_{i+1}^T z_{i+1} / r_{0}^T z_{0}
217  normReduction = sqrt(rdots/dotRr0);
218 
219  Log::info() << "DRPCG end of iteration " << jiter+1 << std::endl;
220  printNormReduction(jiter+1, sqrt(rdots), normReduction);
221  printQuadraticCostFunction(jiter+1, costJ, costJb, costJoJc);
222 
223  // Save the pairs for preconditioning
224  lmp_.push(pp, hh, qq, rho);
225 
226  if (normReduction < tolerance) {
227  Log::info() << "DRPCG: Achieved required reduction in residual norm." << std::endl;
228  break;
229  }
230 
231  vvecs.push_back(rr);
232  zvecs.push_back(zz);
233  scals.push_back(1.0/rdots);
234  }
235 
236 // Generate the (second-level) Limited Memory Preconditioner
237  lmp_.update(B);
238 
239  return normReduction;
240 }
241 
242 // -----------------------------------------------------------------------------
243 
244 } // namespace oops
245 
246 #endif // OOPS_ASSIMILATION_DRPCGMINIMIZER_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
DR (Derber and Rosati) Minimizers.
Definition: DRMinimizer.h:46
DRPCG Minimizer.
QNewtonLMP< CtrlInc_, Bmat_ > lmp_
CostFunction< MODEL, OBS > CostFct_
ControlIncrement< MODEL, OBS > CtrlInc_
BMatrix< MODEL, OBS > Bmat_
HtRinvHMatrix< MODEL, OBS > HtRinvH_
const std::string classname() const override
DRPCGMinimizer(const eckit::Configuration &, const CostFct_ &)
double solve(CtrlInc_ &, CtrlInc_ &, CtrlInc_ &, const Bmat_ &, const HtRinvH_ &, const double, const double, const int, const double) override
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)