OOPS
RPLanczosMinimizer.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_RPLANCZOSMINIMIZER_H_
12 #define OOPS_ASSIMILATION_RPLANCZOSMINIMIZER_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 
29 namespace oops {
30 
31 /// RLanczos Minimizer
32 /*!
33  * \brief Augmented Restricted Lanczos method. It is the Lanczos version
34  * of RPCG.
35  *
36  * This solver is based on the algorithm from Gurol, PhD Manuscript,
37  * 2013. It performs minimization in the observation space.
38  *
39  * It solves the linear system \f$ (I + R^{-1}HBH^T) lambda = H^T R^{-1}d \f$
40  * with \f$ HBH^T \f$ inner-product in the augmented observation space.
41  *
42  * Note that the traditional \f$ B\f$-preconditioning in model space
43  * corresponds to \f$I\f$ for this algorithm.
44  *
45  * A second-level preconditioner, \f$ G \f$, must be symmetric and
46  * positive definite with respect to \f$ HBH^T \f$ inner-product.
47  * Possible preconditioning is detailed in S. Gurol, PhD Manuscript, 2013.
48  *
49  * On entry:
50  * - vv = starting point, \f$ v_0 \f$.
51  * - vvp = starting point augmented part, \f$ vp_0 \f$
52  * - rr = right hand side.
53  * - aa = \f$ H (xb - xk) \f$
54  * - sigma = \f$ (xb - xk)^T B^{-1} (xb - xk)\f$.
55  * - HBHt = \f$ HBH^T \f$.
56  * - Rinv = \f$ R^{-1} \f$.
57  * - G = preconditioner \f$ G \f$.
58  *
59  * On exit, vv and vvp will contain the solution \f$ lambda = [vv; vp] \f$
60  *
61  * The return value is the achieved reduction in preconditioned residual norm.
62  *
63  * Iteration will stop if the maximum iteration limit "maxiter" is reached
64  * or if the preconditioned residual norm reduces by a factor of "tolerance".
65  */
66 
67 // -----------------------------------------------------------------------------
68 
69 template<typename MODEL, typename OBS> class RPLanczosMinimizer : public DualMinimizer<MODEL, OBS> {
74 
75  public:
76  const std::string classname() const override {return "RPLanczosMinimizer";}
77  RPLanczosMinimizer(const eckit::Configuration &, const CostFct_ & J)
78  : DualMinimizer<MODEL, OBS>(J) {}
80 
81  private:
82  double solve(Dual_ &, double &, Dual_ &, const HBHt_ &, const Rinv_ &,
83  const int &, const double &, Dual_ &, const double &) override;
84 };
85 
86 // =============================================================================
87 
88 template<typename MODEL, typename OBS>
89 double RPLanczosMinimizer<MODEL, OBS>::solve(Dual_ & vv, double & vvp, Dual_ & rr,
90  const HBHt_ & HBHt, const Rinv_ & Rinv,
91  const int & maxiter, const double & tolerance,
92  Dual_ & dy, const double & sigma) {
93  IdentityMatrix<Dual_> precond;
94 
95  Dual_ zz(vv);
96  Dual_ ww(vv);
97  Dual_ tt(vv);
98  Dual_ vold(vv);
99  Dual_ v(vv);
100 
101  // augmented part of the vectors
102  double rrp = 1.0;
103  double zzp;
104  double wwp;
105  double ttp = 0.0;
106  double voldp = 0.0;
107  double vp = 1.0;
108 
109  std::vector<Dual_> vVEC; // required for re-orthogonalization
110  std::vector<Dual_> tVEC; // required for re-orthogonalization
111  std::vector<Dual_> zVEC; // required for solution
112  // reserve space to avoid extra copies
113  vVEC.reserve(maxiter+1);
114  tVEC.reserve(maxiter+1);
115  zVEC.reserve(maxiter+1);
116 
117  std::vector<double> vpVEC;
118  std::vector<double> tpVEC;
119  std::vector<double> zpVEC;
120  std::vector<double> alphas;
121  std::vector<double> betas;
122  std::vector<double> y;
123  std::vector<double> d;
124 
125  // zzaug = Gaug rraug
126  precond.multiply(rr, zz);
127  zzp = rrp;
128 
129  // ttaug = HBHtaug zzaug
130  HBHt.multiply(zz, tt);
131  tt.axpy(zzp, dy);
132  ttp = dot_product(dy, zz) + sigma*zzp;
133 
134  double normReduction = 1.0;
135  double beta0 = sqrt(dot_product(rr, tt) + rrp*ttp);
136  double beta = 0.0;
137 
138  vold.zero();
139  v = rr;
140  vp = rrp;
141  v *= 1/beta0;
142  vp *= 1/beta0;
143  tt *= 1/beta0;
144  zz *= 1/beta0;
145  ttp *= 1/beta0;
146  zzp *= 1/beta0;
147 
148  vVEC.clear();
149  zVEC.clear();
150  tVEC.clear();
151  vpVEC.clear();
152  zpVEC.clear();
153  tpVEC.clear();
154 
155  betas.clear();
156  alphas.clear();
157 
158  vVEC.push_back(v);
159  zVEC.push_back(zz);
160  tVEC.push_back(tt);
161  vpVEC.push_back(vp);
162  zpVEC.push_back(zzp);
163  tpVEC.push_back(ttp);
164 
165  int jiter = 0;
166  Log::info() << std::endl;
167  while (jiter < maxiter) {
168  Log::info() << "RPLanczos Starting Iteration " << jiter+1 << std::endl;
169 
170  // ww = (RinvHBHt + I) zz - beta * vold
171  Rinv.multiply(tt, ww);
172  ww += zz;
173  wwp = zzp;
174  ww.axpy(-beta, vold);
175  wwp -= beta*voldp;
176 
177  double alpha = dot_product(tt, ww) + ttp*wwp;
178 
179  ww.axpy(-alpha, v); // w = w - alpha * v
180  wwp -= alpha*vp;
181 
182  // Re-orthogonalization
183  for (int iiter = 0; iiter < jiter; ++iiter) {
184  double proj = dot_product(ww, tVEC[iiter]) + wwp * tpVEC[iiter];
185  ww.axpy(-proj, vVEC[iiter]);
186  wwp -= proj * vpVEC[iiter];
187  }
188 
189  // wwaug = Gaug zzaug
190  precond.multiply(ww, zz);
191  zzp = wwp;
192 
193  // ttaug = HBHtaug zzaug
194  HBHt.multiply(zz, tt);
195  tt.axpy(zzp, dy);
196  ttp = dot_product(dy, zz) + sigma*zzp;
197 
198  beta = sqrt(dot_product(tt, ww) + ttp*wwp);
199 
200  vold = v;
201  voldp = vp;
202  v = ww;
203  vp = wwp;
204  v *= 1/beta;
205  vp *= 1/beta;
206  tt *= 1/beta;
207  zz *= 1/beta;
208  ttp *= 1/beta;
209  zzp *= 1/beta;
210 
211  vVEC.push_back(v);
212  zVEC.push_back(zz);
213  tVEC.push_back(tt);
214  vpVEC.push_back(vp);
215  zpVEC.push_back(zzp);
216  tpVEC.push_back(ttp);
217 
218  alphas.push_back(alpha);
219 
220  y.clear();
221  if (jiter == 0) {
222  y.push_back(beta0/alpha);
223  } else {
224  // Solve the tridiagonal system T_jiter y_jiter = beta0 * e_1
225  d.clear();
226  for (int iiter = 0; iiter <= jiter; ++iiter) {
227  d.push_back(beta0*(dot_product(tVEC[0], vVEC[iiter]) + tpVEC[0]*vpVEC[iiter]));
228  }
229  TriDiagSolve(alphas, betas, d, y);
230  }
231 
232  // Gradient norm in precond metric --> sqrt(r't) --> beta * y(jiter)
233  double rznorm = beta*std::abs(y[jiter]);
234 
235  normReduction = rznorm/beta0;
236 
237  betas.push_back(beta);
238 
239  Log::info() << "RPLanczos end of iteration " << jiter+1 << std::endl;
240  printNormReduction(jiter+1, rznorm, normReduction);
241 
242  ++jiter;
243 
244  if (normReduction < tolerance) {
245  Log::info() << "RPLanczos: Achieved required reduction in residual norm." << std::endl;
246  break;
247  }
248  }
249 
250  // Calculate the solution (xh = Binv x)
251  for (int iiter = 0; iiter < jiter; ++iiter) {
252  vv.axpy(y[iiter], zVEC[iiter]);
253  vvp += y[iiter]*zpVEC[iiter];
254  }
255 
256  Log::info() << "RPLanczos: end" << std::endl;
257 
258  return normReduction;
259 }
260 
261 // -----------------------------------------------------------------------------
262 
263 } // namespace oops
264 
265 #endif // OOPS_ASSIMILATION_RPLANCZOSMINIMIZER_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
double solve(Dual_ &, double &, Dual_ &, const HBHt_ &, const Rinv_ &, const int &, const double &, Dual_ &, const double &) override
RPLanczosMinimizer(const eckit::Configuration &, const CostFct_ &J)
RinvMatrix< MODEL, OBS > Rinv_
CostFunction< MODEL, OBS > CostFct_
const std::string classname() const override
HBHtMatrix< MODEL, OBS > HBHt_
DualVector< MODEL, OBS > Dual_
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 TriDiagSolve(const std::vector< double > &diag, const std::vector< double > &sub, const std::vector< double > &rhs, std::vector< double > &sol)
void printNormReduction(int iteration, const double &grad, const double &norm)