OOPS
RPCGMinimizer.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_RPCGMINIMIZER_H_
12 #define OOPS_ASSIMILATION_RPCGMINIMIZER_H_
13 
14 #include <string>
15 #include <vector>
16 
24 #include "oops/util/dot_product.h"
25 #include "oops/util/formats.h"
26 #include "oops/util/Logger.h"
27 
28 namespace oops {
29 
30 /// RPCG Minimizer
31 /*!
32  * \brief Augmented Restricted Preconditioned Conjugate Gradients.
33  *
34  * This solver is based on the algorithm proposed in Gratton and Tshimanga,
35  * QJRMS, 135: 1573-1585 (2009). It performs minimization in observation space.
36  *
37  * It solves the linear system \f$ (I + R^{-1}HBH^T) lambda = H^T R^{-1}d \f$
38  * with \f$ HBH^T \f$ inner-product in the augmented observation space.
39  *
40  * Note that the traditional \f$ B\f$-preconditioning in model space
41  * corresponds to \f$I\f$ for this algorithm.
42  *
43  * A second-level preconditioner, \f$ G \f$, must be symmetric and
44  * positive definite with respect to \f$ HBH^T \f$ inner-product.
45  * Possible preconditioning is detailed in S. Gurol, PhD Manuscript, 2013.
46  *
47  * On entry:
48  * - vv = starting point, \f$ v_0 \f$.
49  * - vvp = starting point augmented part, \f$ vp_0 \f$
50  * - rr = right hand side.
51  * - dy = \f$ H (xb - xk) \f$
52  * - sigma = \f$ (xb - xk)^T B^{-1} (xb - xk)\f$.
53  * - HBHt = \f$ HBH^T \f$.
54  * - Rinv = \f$ R^{-1} \f$.
55  * - G = preconditioner \f$ G \f$.
56  * - Gt = preconditioner transpose \f$ G^T \f$.
57  *
58  * On exit, vv and vvp will contain the solution \f$ lambda = [vv; vp] \f$
59  *
60  * The return value is the achieved reduction in preconditioned residual norm.
61  *
62  * Iteration will stop if the maximum iteration limit "maxiter" is reached
63  * or if the preconditioned residual norm reduces by a factor of "tolerance".
64  */
65 
66 // -----------------------------------------------------------------------------
67 
68 template<typename MODEL, typename OBS> class RPCGMinimizer : public DualMinimizer<MODEL, OBS> {
73 
74  public:
75  const std::string classname() const override {return "RPCGMinimizer";}
76  RPCGMinimizer(const eckit::Configuration &, const CostFct_ & J): DualMinimizer<MODEL, OBS>(J) {}
78 
79  private:
80  double solve(Dual_ &, double &, Dual_ &, const HBHt_ &, const Rinv_ &,
81  const int &, const double &, Dual_ &, const double &) override;
82 };
83 
84 // =============================================================================
85 
86 template<typename MODEL, typename OBS>
87 double RPCGMinimizer<MODEL, OBS>::solve(Dual_ & vv, double & vvp, Dual_ & rr,
88  const HBHt_ & HBHt, const Rinv_ & Rinv,
89  const int & maxiter, const double & tolerance,
90  Dual_ & dy, const double & sigma) {
91  IdentityMatrix<Dual_> precond;
92  IdentityMatrix<Dual_> precondt;
93 
94  Dual_ zz(vv);
95  Dual_ ww(vv);
96  Dual_ pp(vv);
97  Dual_ tt(vv);
98  Dual_ ll(vv);
99  Dual_ qq(vv);
100  Dual_ w(vv);
101  Dual_ v(vv);
102 
103  // augmented part of the vectors
104  double rrp = 1.0;
105  double llp;
106  double wwp;
107  double zzp;
108  double ppp = 0.0;
109  double ttp = 0.0;
110  double qqp;
111  double vp = 1.0;
112  double wp;
113 
114  std::vector<Dual_> vVEC; // required for re-orthogonalization
115  std::vector<Dual_> wVEC; // required for re-orthogonalization
116  // reserve space to avoid extra copies
117  vVEC.reserve(maxiter+1);
118  wVEC.reserve(maxiter+1);
119 
120  std::vector<double> vpVEC;
121  std::vector<double> wpVEC;
122 
123  // llaug = HBHtaug rraug
124  HBHt.multiply(rr, ll);
125  ll.axpy(rrp, dy);
126  llp = dot_product(dy, rr) + sigma*rrp;
127 
128  // wwaug = Gtaug llaug
129  precondt.multiply(ll, ww); // UPDATE WHEN G IS NOT IDENTITY
130  wwp = llp;
131 
132  // zzaug = Gaug rraug
133  precond.multiply(rr, zz); // UPDATE WHEN G IS NOT IDENTITY
134  zzp = rrp;
135 
136  double dotwr = dot_product(rr, ww) + rrp*wwp;
137  double normReduction = 1.0;
138  double dotw0r0 = dotwr;
139  double dotwr_old = 0.0;
140 
141  v = rr;
142  v *= 1/sqrt(dotwr);
143  vp = rrp;
144  vp *= 1/sqrt(dotwr);
145  w = ww;
146  w *= 1/sqrt(dotwr);
147  wp = wwp;
148  wp *= 1/sqrt(dotwr);
149 
150  vVEC.clear();
151  vpVEC.clear();
152  wVEC.clear();
153  wpVEC.clear();
154 
155  vVEC.push_back(v);
156  vpVEC.push_back(vp);
157  wVEC.push_back(w);
158  wpVEC.push_back(wp);
159 
160  Log::info() << std::endl;
161  for (int jiter = 0; jiter < maxiter; ++jiter) {
162  Log::info() << "RPCG Starting Iteration " << jiter+1 << std::endl;
163 
164  if (jiter > 0) {
165  double beta = dotwr/dotwr_old;
166  Log::info() << "RPCG beta = " << beta << std::endl;
167 
168  pp *= beta;
169  ppp *= beta;
170 
171  tt *= beta;
172  ttp *= beta;
173  }
174 
175  pp += zz;
176  ppp += zzp; // ppaug = zzaug + beta*ppaug
177 
178  tt += ww;
179  ttp += wwp; // ttaug = wwaug + beta*ttaug
180 
181  // (RinvHBHt + I) pp
182  Rinv.multiply(tt, qq);
183  qq += pp;
184  qqp = ppp; // qqaug = ppaug + Rinv_aug ttaug
185 
186  double alpha = dotwr/(dot_product(qq, tt) + qqp * ttp);
187  Log::info() << "RPCG alpha = " << alpha << std::endl;
188 
189  vv.axpy(alpha, pp);
190  vvp = vvp + alpha*ppp; // vvaug = vvaug + alpha*ppaug
191 
192  rr.axpy(-alpha, qq);
193  rrp = rrp - alpha*qqp; // rraug = rraug - alpha*qqaug
194 
195  // Re-orthogonalization
196  for (int iiter = 0; iiter < jiter; ++iiter) {
197  double proj = dot_product(rr, wVEC[iiter]) + rrp * wpVEC[iiter];
198  rr.axpy(-proj, vVEC[iiter]);
199  rrp -= proj*vpVEC[iiter];
200  }
201 
202  // llaug = HBHt_aug rraug
203  HBHt.multiply(rr, ll);
204  ll.axpy(rrp, dy);
205  llp = dot_product(dy, rr) + sigma * rrp;
206 
207  // wwaug = Gt_aug llaug
208  precondt.multiply(ll, ww); // UPDATE WHEN G IS NOT IDENTITY
209  wwp = llp;
210 
211  // zzaug = Gaug rraug
212  precond.multiply(rr, zz); // UPDATE WHEN G IS NOT IDENTITY
213  zzp = rrp;
214 
215  dotwr_old = dotwr;
216  dotwr = dot_product(ww, rr) + rrp * wwp;
217 
218  v = rr;
219  v *= 1/sqrt(dotwr);
220  vp = rrp;
221  vp *= 1/sqrt(dotwr);
222  w = ww;
223  w *= 1/sqrt(dotwr);
224  wp = wwp;
225  wp *= 1/sqrt(dotwr);
226 
227  vVEC.push_back(v);
228  vpVEC.push_back(vp);
229  wVEC.push_back(w);
230  wpVEC.push_back(wp);
231 
232  normReduction = sqrt(dotwr/dotw0r0);
233 
234  Log::info() << "RPCG end of iteration " << jiter+1 << std::endl;
235  printNormReduction(jiter+1, sqrt(dotwr), normReduction);
236 
237  if (normReduction < tolerance) {
238  Log::info() << "RPCG: Achieved required reduction in residual norm." << std::endl;
239  break;
240  }
241  }
242  return normReduction;
243 }
244 
245 // -----------------------------------------------------------------------------
246 
247 } // namespace oops
248 
249 #endif // OOPS_ASSIMILATION_RPCGMINIMIZER_H_
Cost Function.
Definition: CostFunction.h:53
Dual Minimizer.
Definition: DualMinimizer.h:40
Container of dual space vectors for all terms of the cost function.
Definition: DualVector.h:35
void axpy(const double, const DualVector &)
Definition: DualVector.h:215
The matrix.
Definition: HBHtMatrix.h:35
void multiply(const Dual_ &dy, Dual_ &dz) const
Definition: HBHtMatrix.h:61
Identity matrix.
void multiply(const VECTOR &a, VECTOR &b) const
RPCG Minimizer.
Definition: RPCGMinimizer.h:68
RPCGMinimizer(const eckit::Configuration &, const CostFct_ &J)
Definition: RPCGMinimizer.h:76
DualVector< MODEL, OBS > Dual_
Definition: RPCGMinimizer.h:70
HBHtMatrix< MODEL, OBS > HBHt_
Definition: RPCGMinimizer.h:71
CostFunction< MODEL, OBS > CostFct_
Definition: RPCGMinimizer.h:69
RinvMatrix< MODEL, OBS > Rinv_
Definition: RPCGMinimizer.h:72
double solve(Dual_ &, double &, Dual_ &, const HBHt_ &, const Rinv_ &, const int &, const double &, Dual_ &, const double &) override
Definition: RPCGMinimizer.h:87
const std::string classname() const override
Definition: RPCGMinimizer.h:75
The matrix.
Definition: RinvMatrix.h:29
void multiply(const Dual_ &dx, Dual_ &dy) const
Definition: RinvMatrix.h:36
The namespace for the main oops code.
void printNormReduction(int iteration, const double &grad, const double &norm)