OOPS
DRMinimizer.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_DRMINIMIZER_H_
12 #define OOPS_ASSIMILATION_DRMINIMIZER_H_
13 
14 #include <limits>
15 #include <memory>
16 #include <string>
17 #include <vector>
18 
19 #include "eckit/config/Configuration.h"
26 
27 #include "oops/util/dot_product.h"
28 #include "oops/util/FloatCompare.h"
29 #include "oops/util/formats.h"
30 #include "oops/util/Logger.h"
31 
32 namespace oops {
33 
34 /// DR (Derber and Rosati) Minimizers
35 /*!
36  * DRMinimizer is the base class for all minimizers that use \f$ B\f$ to
37  * precondition the variational minimisation problem and use the auxiliary
38  * variable \f$ \hat{x}=B^{-1}x\f$ and to update it in parallel to \f$ x\f$
39  * based on Derber and Rosati, 1989, J. Phys. Oceanog. 1333-1347.
40  * \f$ J_b\f$ is then computed as \f$ x^T\hat{x}\f$ eliminating the need for
41  * \f$ B^{-1}\f$ or \f$ B^{1/2}\f$.
42  */
43 
44 // -----------------------------------------------------------------------------
45 
46 template<typename MODEL, typename OBS> class DRMinimizer : public Minimizer<MODEL, OBS> {
53 
54  public:
55  explicit DRMinimizer(const CostFct_ & J): Minimizer_(J), J_(J), gradJb_(), costJ0Jb_(0) {}
57  const std::string classname() const override = 0;
58 
59  private:
60  CtrlInc_ * doMinimize(const eckit::Configuration &) override;
61  virtual double solve(CtrlInc_ &, CtrlInc_ &, CtrlInc_ &,
62  const Bmat_ &, const HtRinvH_ &,
63  const double, const double,
64  const int, const double) = 0;
65 
66  const CostFct_ & J_;
67  std::unique_ptr<CtrlInc_> gradJb_;
68  std::vector<CtrlInc_> dxh_;
69  double costJ0Jb_;
70 };
71 
72 // =============================================================================
73 
74 template<typename MODEL, typename OBS>
76 DRMinimizer<MODEL, OBS>::doMinimize(const eckit::Configuration & config) {
77  int ninner = config.getInt("ninner");
78  double gnreduc = config.getDouble("gradient norm reduction");
79 
80  bool runOnlineAdjTest = config.getBool("online diagnostics.online adj test", false);
81 
82  if (gradJb_) {
83  gradJb_.reset(new CtrlInc_(J_.jb().resolution(), *gradJb_));
84  } else {
85  gradJb_.reset(new CtrlInc_(J_.jb()));
86  }
87 
88  Log::info() << std::endl;
89  Log::info() << classname() << ": max iter = " << ninner
90  << ", requested norm reduction = " << gnreduc << std::endl;
91 
92 // Define the matrices
93  const Bmat_ B(J_);
94  const HtRinvH_ HtRinvH(J_, runOnlineAdjTest);
95 
96 // Define minimisation starting point
97  // dx
98  CtrlInc_ * dx = new CtrlInc_(J_.jb());
99  // dxh = B^{-1} dx
100  CtrlInc_ dxh(J_.jb());
101 
102 // Compute RHS (sum B^{-1} dx_{i}) + H^T R^{-1} d
103 // dx_i = x_i - x_{i-1}; dx_1 = x_1 - x_b
104  CtrlInc_ rhs(J_.jb());
105  CtrlInc_ dx0(rhs);
106  if (config.has("fsoi")) {
107  const eckit::LocalConfiguration FcSensitivityConfig(config, "fsoi.input forecast sensitivity");
108  rhs.read(FcSensitivityConfig);
109  dx0 = rhs;
110  dx0 *= -1.0;
111  Log::info() << classname() << " rhs has forecast sensitivity" << std::endl;
112  } else {
113  J_.computeGradientFG(rhs);
114  J_.jb().addGradientFG(rhs, *gradJb_);
115  }
116  rhs *= -1.0;
117  Log::info() << classname() << " rhs" << rhs << std::endl;
118 
119 // Check for zero gradient (for example if no obs)
120  const double gnorm = dot_product(rhs, rhs);
121  const double epsilon = config.getDouble("epsilon", std::numeric_limits<double>::epsilon());
122  Log::info() << "Initial RHS squared norm = " << gnorm << std::endl;
123  if (gnorm < epsilon) {
124  Log::info() << "RHS smaller than " << epsilon << ", returning." << std::endl;
125  return dx;
126  }
127 
128 // Set J[0] = 0.5 (x_i - x_b)^T B^{-1} (x_i - x_b) + 0.5 d^T R^{-1} d
129  const double costJ0Jb = costJ0Jb_;
130  const double costJ0JoJc = J_.getCostJoJc();
131 
132 // Solve the linear system
133  double reduc = this->solve(*dx, dxh, rhs, B, HtRinvH, costJ0Jb, costJ0JoJc, ninner, gnreduc);
134 
135  Log::test() << classname() << ": reduction in residual norm = " << reduc << std::endl;
136  Log::info() << classname() << " output increment:" << *dx << std::endl;
137 
138 // Update gradient Jb
139  *gradJb_ += dxh;
140  dxh_.push_back(dxh);
141 
142 // Update Jb component of J[0]: 0.5 (x_i - x_b)^T B^-1 (x_i - x_b)
143  costJ0Jb_ += 0.5 * dot_product(*dx, dxh);
144  for (unsigned int jouter = 1; jouter < dxh_.size(); ++jouter) {
145  CtrlInc_ dxhtmp(dx->geometry(), dxh_[jouter-1]);
146  costJ0Jb_ += dot_product(*dx, dxhtmp);
147  }
148 
149  if (config.has("fsoi")) {
150  Log::info() << classname() << " Entering Observation Sensitivity Calculation" << std::endl;
151 
152  // Multiply result of solver by RinvH to get observation sensitivity (ys)
154  const RinvH_ RinvH(J_);
155  RinvH.multiply(*dx, ys);
156 
157  // Write out observation sensitivity
158  const std::string osensname = "ObsSensitivity";
159  ys.saveDep(osensname);
160 
161  bool runFSOIincTest = config.getBool("fsoi.increment test", false);
162  if (runFSOIincTest) {
163  // Get departures
165  for (unsigned jj = 0; jj < J_.nterms(); ++jj) {
166  std::unique_ptr<GeneralizedDepartures> ww(J_.jterm(jj).newGradientFG());
167  dp.append(J_.jterm(jj).multiplyCovar(*ww));
168  }
169 
170  // <K dp, dx>, where dx = K dp
171  double adj_tst_fwd = dot_product(dx0, dx0);
172  // <dp, Kt dx>, where K = Hessian Ht Rinv; dp=departures
173  double adj_tst_bwd = dot_product(ys, dp);
174 
175  Log::info() << "Online FSOI increment test: " << std::endl
176  << util::PrintAdjTest(adj_tst_fwd, adj_tst_bwd, "K") << std::endl;
177 
178  double fsoi_inctest_tolerance = config.getDouble("fsoi.increment test tolerance", 1.0e-5);
179  bool passed = oops::is_close_absolute(adj_tst_fwd, adj_tst_bwd, fsoi_inctest_tolerance);
180  if (passed) {
181  Log::test() << "FSOI increment test within tolerance." << std::endl;
182  } else {
183  Log::test() << "FSOI increment test fails tolerance bound." << std::endl;
184  }
185  }
186 
187  // Make sure not to update state in FSOI mode
188  dx->zero();
189  }
190 
191  return dx;
192 }
193 
194 // -----------------------------------------------------------------------------
195 
196 } // namespace oops
197 
198 #endif // OOPS_ASSIMILATION_DRMINIMIZER_H_
The matrix.
Definition: BMatrix.h:27
void zero()
Linear algebra operators.
Geometry_ geometry() const
Get geometry.
void read(const eckit::Configuration &)
I/O and diagnostics.
Cost Function.
Definition: CostFunction.h:53
DR (Derber and Rosati) Minimizers.
Definition: DRMinimizer.h:46
std::unique_ptr< CtrlInc_ > gradJb_
Definition: DRMinimizer.h:67
CostFunction< MODEL, OBS > CostFct_
Definition: DRMinimizer.h:48
HtRinvHMatrix< MODEL, OBS > HtRinvH_
Definition: DRMinimizer.h:50
BMatrix< MODEL, OBS > Bmat_
Definition: DRMinimizer.h:47
const CostFct_ & J_
Definition: DRMinimizer.h:66
ControlIncrement< MODEL, OBS > CtrlInc_
Definition: DRMinimizer.h:49
CtrlInc_ * doMinimize(const eckit::Configuration &) override
Definition: DRMinimizer.h:76
const std::string classname() const override=0
Minimizer< MODEL, OBS > Minimizer_
Definition: DRMinimizer.h:51
std::vector< CtrlInc_ > dxh_
Definition: DRMinimizer.h:68
virtual double solve(CtrlInc_ &, CtrlInc_ &, CtrlInc_ &, const Bmat_ &, const HtRinvH_ &, const double, const double, const int, const double)=0
RinvHMatrix< MODEL, OBS > RinvH_
Definition: DRMinimizer.h:52
DRMinimizer(const CostFct_ &J)
Definition: DRMinimizer.h:55
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
void saveDep(const std::string &) const
Definition: DualVector.h:253
A Minimizer knows how to minimize a cost function.
Definition: Minimizer.h:37
The matrix.
Definition: RinvHMatrix.h:38
void multiply(const CtrlInc_ &dx, Dual_ &zz) const
Definition: RinvHMatrix.h:62
The namespace for the main oops code.