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