OOPS
DualMinimizer.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_DUALMINIMIZER_H_
12 #define OOPS_ASSIMILATION_DUALMINIMIZER_H_
13 
14 #include <limits>
15 #include <memory>
16 #include <string>
17 
18 #include "eckit/config/Configuration.h"
28 #include "oops/util/dot_product.h"
29 #include "oops/util/Logger.h"
30 
31 namespace oops {
32 
33 /// Dual Minimizer
34 /*!
35  * Base class for all dual (observation) space minimizers.
36  */
37 
38 // -----------------------------------------------------------------------------
39 
40 template<typename MODEL, typename OBS> class DualMinimizer : public Minimizer<MODEL, OBS> {
48 
49  public:
50  explicit DualMinimizer(const CostFct_ & J): Minimizer_(J), J_(J), gradJb_() {}
52  const std::string classname() const override = 0;
53 
54  private:
55  CtrlInc_ * doMinimize(const eckit::Configuration &) override;
56  virtual double solve(Dual_ &, double &, Dual_ &, const HBHt_ &, const Rinv_ &,
57  const int &, const double &, Dual_ &, const double &) = 0;
58 
59  const CostFct_ & J_;
60  std::unique_ptr<CtrlInc_> gradJb_;
61 };
62 
63 // =============================================================================
64 
65 template<typename MODEL, typename OBS>
67 DualMinimizer<MODEL, OBS>::doMinimize(const eckit::Configuration & config) {
68  int ninner = config.getInt("ninner");
69  double gnreduc = config.getDouble("gradient norm reduction");
70 
71  bool runOnlineAdjTest = config.getBool("online diagnostics.online adj test", false);
72 
73  if (gradJb_) {
74  gradJb_.reset(new CtrlInc_(J_.jb().resolution(), *gradJb_));
75  } else {
76  gradJb_.reset(new CtrlInc_(J_.jb()));
77  }
78 
79  Log::info() << std::endl;
80  Log::info() << classname() << ": max iter = " << ninner
81  << ", requested norm reduction = " << gnreduc << std::endl;
82 
83 // Define the matrices
84  const Bmat_ B(J_);
85  const HBHt_ HBHt(J_, runOnlineAdjTest);
86  const Rinv_ Rinv(J_);
87  const HMatrix<MODEL, OBS> H(J_);
88  const HtMatrix<MODEL, OBS> Ht(J_);
89 
90  CtrlInc_ * dx = new CtrlInc_(J_.jb());
91 
92 // Define minimisation starting point in dual space
93  Dual_ vv;
94  for (unsigned jj = 0; jj < J_.nterms(); ++jj) {
95  vv.append(J_.jterm(jj).newDualVector());
96  }
97  double vvp = 0.0;
98 
99 // Get R^{-1} d
100  Dual_ rr;
101  for (unsigned jj = 0; jj < J_.nterms(); ++jj) {
102  rr.append(J_.jterm(jj).newGradientFG());
103  }
104  rr *= -1.0;
105 
106 // Check for zero gradient (for example if no obs)
107  const double rnorm = dot_product(rr, rr);
108  const double epsilon = config.getDouble("epsilon", std::numeric_limits<double>::epsilon());
109  Log::info() << "Initial RHS squared norm = " << rnorm << std::endl;
110  if (rnorm < epsilon) {
111  Log::info() << "RHS smaller than " << epsilon << ", returning." << std::endl;
112  return dx;
113  }
114 
115 // Update rr if initial dx in model space is not a zero vector
116 // rr = rr - Rinv H dx0
117 
118 // Compute a = H (xb - xk)
119  Dual_ dy;
120  CtrlInc_ hdxfg(J_.jb().getFirstGuess());
121  H.multiply(hdxfg, dy);
122  dy *= -1.0;
123 
124 // Compute sigma = (xb - xk)^T Binv (xb - xk)
125  CtrlInc_ g0(J_.jb());
126  J_.jb().addGradientFG(g0, *gradJb_);
127 
128  double sigma = dot_product(J_.jb().getFirstGuess(), g0);
129 
130 // Solve the linear system
131  double reduc = this->solve(vv, vvp, rr, HBHt, Rinv, ninner, gnreduc, dy, sigma);
132 
133  Log::test() << classname() << ": reduction in residual norm = " << reduc << std::endl;
134 
135 // Recover solution in primal space
136  CtrlInc_ dh(J_.jb());
137  J_.zeroAD(dh);
138  Ht.multiply(vv, dh);
139  dh.axpy(-vvp, g0);
140 
141  B.multiply(dh, *dx); // BHtaug vvaug
142 
143  Log::info() << classname() << ": Estimated Final Jb = "
144  << 0.5 * dot_product(*dx, dh) << std::endl;
145  Log::info() << classname() << " output" << *dx << std::endl;
146 
147 // Update gradient Jb
148  *gradJb_ += dh;
149 
150  return dx;
151 }
152 
153 // -----------------------------------------------------------------------------
154 
155 } // namespace oops
156 
157 #endif // OOPS_ASSIMILATION_DUALMINIMIZER_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
Dual Minimizer.
Definition: DualMinimizer.h:40
DualVector< MODEL, OBS > Dual_
Definition: DualMinimizer.h:44
const std::string classname() const override=0
CtrlInc_ * doMinimize(const eckit::Configuration &) override
Definition: DualMinimizer.h:67
RinvMatrix< MODEL, OBS > Rinv_
Definition: DualMinimizer.h:47
DualMinimizer(const CostFct_ &J)
Definition: DualMinimizer.h:50
std::unique_ptr< CtrlInc_ > gradJb_
Definition: DualMinimizer.h:60
Minimizer< MODEL, OBS > Minimizer_
Definition: DualMinimizer.h:46
CostFunction< MODEL, OBS > CostFct_
Definition: DualMinimizer.h:42
HBHtMatrix< MODEL, OBS > HBHt_
Definition: DualMinimizer.h:45
const CostFct_ & J_
Definition: DualMinimizer.h:59
ControlIncrement< MODEL, OBS > CtrlInc_
Definition: DualMinimizer.h:41
BMatrix< MODEL, OBS > Bmat_
Definition: DualMinimizer.h:43
virtual double solve(Dual_ &, double &, Dual_ &, const HBHt_ &, const Rinv_ &, const int &, const double &, Dual_ &, const double &)=0
Container of dual space vectors for all terms of the cost function.
Definition: DualVector.h:35
void append(std::unique_ptr< GeneralizedDepartures > &&)
Definition: DualVector.h:110
The matrix.
Definition: HBHtMatrix.h:35
The matrix.
Definition: HMatrix.h:36
void multiply(CtrlInc_ &dx, DualVector< MODEL, OBS > &dy, const bool idModel=false) const
Definition: HMatrix.h:44
The matrix.
Definition: HtMatrix.h:33
void multiply(const DualVector< MODEL, OBS > &dy, ControlIncrement< MODEL, OBS > &dx, const bool idModel=false) const
Definition: HtMatrix.h:40
A Minimizer knows how to minimize a cost function.
Definition: Minimizer.h:37
The matrix.
Definition: RinvMatrix.h:29
The namespace for the main oops code.