OOPS
PrimalMinimizer.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_PRIMALMINIMIZER_H_
12 #define OOPS_ASSIMILATION_PRIMALMINIMIZER_H_
13 
14 #include <limits>
15 #include <memory>
16 #include <string>
17 
18 #include "eckit/config/Configuration.h"
24 #include "oops/util/FloatCompare.h"
25 #include "oops/util/Logger.h"
26 
28 
29 namespace oops {
30 
31 /// Primal Minimizer
32 /*!
33  * PrimalMinimizer is the base class for all minimizers that minimize the
34  * variational data assimilation cost function in primal (model) space.
35  */
36 
37 // -----------------------------------------------------------------------------
38 
39 template<typename MODEL, typename OBS> class PrimalMinimizer : public Minimizer<MODEL, OBS> {
46 
47  public:
48  explicit PrimalMinimizer(const CostFct_ & J): Minimizer_(J), J_(J) {}
50  const std::string classname() const override = 0;
51 
52  private:
53  CtrlInc_ * doMinimize(const eckit::Configuration &) override;
54  virtual double solve(CtrlInc_ &, const CtrlInc_ &,
55  const Hessian_ &, const Bmat_ &,
56  const int, const double) = 0;
57 
58  const CostFct_ & J_;
59 };
60 
61 // =============================================================================
62 
63 template<typename MODEL, typename OBS>
65 PrimalMinimizer<MODEL, OBS>::doMinimize(const eckit::Configuration & config) {
66  int ninner = config.getInt("ninner");
67  double gnreduc = config.getDouble("gradient norm reduction");
68 
69  bool runOnlineAdjTest = config.getBool("online diagnostics.online adj test", false);
70 
71  Log::info() << classname() << ": max iter = " << ninner
72  << ", requested norm reduction = " << gnreduc << std::endl;
73 
74 // Define the matrices
75  Hessian_ hessian(J_, runOnlineAdjTest);
76  Bmat_ B(J_);
77 
78 // Define minimisation starting point
79  CtrlInc_ * dx = new CtrlInc_(J_.jb());
80 
81 // Compute RHS
82  CtrlInc_ rhs(J_.jb());
83  if (config.has("fsoi")) {
84  const eckit::LocalConfiguration FcSensitivityConfig(config, "fsoi.input forecast sensitivity");
85  rhs.read(FcSensitivityConfig);
86  Log::info() << classname() << " rhs has forecast sensitivity" << std::endl;
87  } else {
88  J_.computeGradientFG(rhs);
89  J_.jb().addGradientFG(rhs);
90  }
91  rhs *= -1.0;
92  Log::info() << classname() << " rhs" << rhs << std::endl;
93 
94 // Check for zero gradient (for example if no obs)
95  const double gnorm = dot_product(rhs, rhs);
96  const double epsilon = config.getDouble("epsilon", std::numeric_limits<double>::epsilon());
97  Log::info() << "Initial RHS squared norm = " << gnorm << std::endl;
98  if (gnorm < epsilon) {
99  Log::info() << "RHS smaller than " << epsilon << ", returning." << std::endl;
100  return dx;
101  }
102 
103 // Solve the linear system
104  double reduc = this->solve(*dx, rhs, hessian, B, ninner, gnreduc);
105 
106  Log::test() << classname() << ": reduction in residual norm = " << reduc << std::endl;
107  Log::info() << classname() << " output" << *dx << std::endl;
108 
109  if (config.has("fsoi")) {
110  Log::info() << classname() << " Entering Observation Sensitivity Calculation" << std::endl;
111 
112  // Multiply result of solver by RinvH to get observation sensitivity (ys)
114  const RinvH_ RinvH(J_);
115  RinvH.multiply(*dx, ys);
116 
117  // Write out observation sensitivity
118  const std::string osensname = "ObsSensitivity";
119  ys.saveDep(osensname);
120 
121  bool runFSOIincTest = config.getBool("fsoi.increment test", false);
122  if (runFSOIincTest) {
123  // Get departures
125  for (unsigned jj = 0; jj < J_.nterms(); ++jj) {
126  std::unique_ptr<GeneralizedDepartures> ww(J_.jterm(jj).newGradientFG());
127  dp.append(J_.jterm(jj).multiplyCovar(*ww));
128  }
129 
130  // <K dp, dx>, where dx = K dp
131  double adj_tst_fwd = dot_product(rhs, rhs);
132  // <dp, Kt dx>, where K = Hessian Ht Rinv; dp=departures
133  double adj_tst_bwd = dot_product(ys, dp);
134 
135  Log::info() << "Online FSOI increment test: " << std::endl
136  << util::PrintAdjTest(adj_tst_fwd, adj_tst_bwd, "K") << std::endl;
137 
138  double fsoi_inctest_tolerance = config.getDouble("fsoi.increment test tolerance", 1.0e-5);
139  bool passed = oops::is_close_absolute(adj_tst_fwd, adj_tst_bwd, fsoi_inctest_tolerance);
140  if (passed) {
141  Log::test() << "FSOI increment test within tolerance." << std::endl;
142  } else {
143  Log::test() << "FSOI increment test fails tolerance bound." << std::endl;
144  }
145  }
146 
147  // Make sure not to update state in FSOI mode
148  dx->zero();
149  }
150 
151  return dx;
152 }
153 
154 // -----------------------------------------------------------------------------
155 
156 } // namespace oops
157 
158 #endif // OOPS_ASSIMILATION_PRIMALMINIMIZER_H_
The matrix.
Definition: BMatrix.h:27
void zero()
Linear algebra operators.
void read(const eckit::Configuration &)
I/O and diagnostics.
Cost Function.
Definition: CostFunction.h:53
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
The Hessian matrix: .
Definition: HessianMatrix.h:34
A Minimizer knows how to minimize a cost function.
Definition: Minimizer.h:37
Primal Minimizer.
virtual double solve(CtrlInc_ &, const CtrlInc_ &, const Hessian_ &, const Bmat_ &, const int, const double)=0
BMatrix< MODEL, OBS > Bmat_
ControlIncrement< MODEL, OBS > CtrlInc_
const std::string classname() const override=0
CtrlInc_ * doMinimize(const eckit::Configuration &) override
Minimizer< MODEL, OBS > Minimizer_
HessianMatrix< MODEL, OBS > Hessian_
const CostFct_ & J_
CostFunction< MODEL, OBS > CostFct_
RinvHMatrix< MODEL, OBS > RinvH_
PrimalMinimizer(const CostFct_ &J)
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.