OOPS
Minimizer.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_MINIMIZER_H_
12 #define OOPS_ASSIMILATION_MINIMIZER_H_
13 
14 #include <map>
15 #include <string>
16 #include <boost/noncopyable.hpp>
17 
18 #include "eckit/config/Configuration.h"
26 #include "oops/base/State.h"
27 #include "oops/util/dot_product.h"
28 #include "oops/util/formats.h"
29 #include "oops/util/Logger.h"
30 #include "oops/util/PrintAdjTest.h"
31 
32 namespace oops {
33 
34 // -----------------------------------------------------------------------------
35 
36 /// A Minimizer knows how to minimize a cost function
37 template<typename MODEL, typename OBS> class Minimizer : private boost::noncopyable {
44 
45  public:
46  explicit Minimizer(const CostFct_ & J): J_(J), outerIteration_(0) {}
47  virtual ~Minimizer() {}
48  ControlIncrement<MODEL, OBS> * minimize(const eckit::Configuration &);
49  virtual const std::string classname() const = 0;
50 
51  private:
52  virtual ControlIncrement<MODEL, OBS> * doMinimize(const eckit::Configuration &) = 0;
53 
54  void adjTests(const eckit::Configuration &);
55  void adjModelTest(const Ht_ &, const H_ &);
56  void adjObsTest(const Ht_ &, const H_ &);
57 
58  void tlmTests(const eckit::Configuration &);
59  void tlmApproxTest(const H_ &);
60  void tlmTaylorTest(const H_ &);
61 
62  void tlmPropagTest(const eckit::Configuration & config, const CtrlInc_ &);
63 
64  const CostFct_ & J_;
66 };
67 
68 // -----------------------------------------------------------------------------
69 
70 template<typename MODEL, typename OBS>
72 Minimizer<MODEL, OBS>::minimize(const eckit::Configuration & config) {
73  // TLM tests
74  this->tlmTests(config);
75 
76  // ADJ tests
77  this->adjTests(config);
78 
79  // Minimize
80  ControlIncrement<MODEL, OBS> * dx = this->doMinimize(config);
81 
82  // Write increment
83  writeIncrement(config, *dx, outerIteration_);
84 
85  // TLM propagation test
86  this->tlmPropagTest(config, *dx);
87 
88  // Update outer loop counter
89  outerIteration_++;
90 
91  return dx;
92 }
93 
94 // -----------------------------------------------------------------------------
95 
96 template<typename MODEL, typename OBS>
97 void Minimizer<MODEL, OBS>::tlmTests(const eckit::Configuration & config) {
98 // Tangent Linear tests
99 
100  if (config.has("online diagnostics")) {
101  const eckit::LocalConfiguration onlineDiag(config, "online diagnostics");
102  bool runTlmTaylorTest = onlineDiag.getBool("tlm taylor test", false);
103  bool runTlmApproxTest = onlineDiag.getBool("tlm approx test", false);
104 
105  const H_ H(J_);
106  const Ht_ Ht(J_);
107 
108  // Online tests
109  // TLM linear approximation test
110  if (outerIteration_ == 0 && runTlmApproxTest) this->tlmApproxTest(H);
111 
112  // TLM Taylor test
113  if (outerIteration_ == 0 && runTlmTaylorTest) this->tlmTaylorTest(H);
114  }
115 }
116 
117 // -----------------------------------------------------------------------------
118 
119 template<typename MODEL, typename OBS>
120 void Minimizer<MODEL, OBS>::adjTests(const eckit::Configuration & config) {
121 // Adjoint tests
122 
123  if (config.has("online diagnostics")) {
124  const eckit::LocalConfiguration onlineDiag(config, "online diagnostics");
125  bool runAdjTlmTest = onlineDiag.getBool("adj tlm test", false);
126  bool runAdjObsTest = onlineDiag.getBool("adj obs test", false);
127 
128  const H_ H(J_);
129  const Ht_ Ht(J_);
130 
131  // Online tests
132  // Model adjoint test
133  if (runAdjTlmTest) this->adjModelTest(Ht, H);
134 
135  // Obs adjoint test
136  if (runAdjObsTest) this->adjObsTest(Ht, H);
137  }
138 }
139 
140 // -----------------------------------------------------------------------------
141 
142 template<typename MODEL, typename OBS>
144 /* TL approx test:
145  calculate and store:
146  M(x + dx), M(x), M'(dx)
147  where dx is initialized using randomization
148 
149  TO DO:
150  write to file depending on requriements */
151 
152  Log::info() << std::endl
153  << "TLM Linear Approximation Test - starting: " << outerIteration_
154  << std::endl << std::endl;
155 
156 // initialize incremets
157  CtrlInc_ dx(J_.jb());
158 
159 // randomize increments
160  J_.jb().randomize(dx);
161 
162 // run NL for unperturbed initial condition
164  ControlVariable<MODEL, OBS> mxx(J_.jb().getBackground());
165  J_.runNL(mxx, pp);
166 
167 // run NL for perturbed initial condition
168  ControlVariable<MODEL, OBS> mpertxx(J_.jb().getBackground());
169  J_.addIncrement(mpertxx, dx);
170  J_.runNL(mpertxx, pp);
171 
172 // run TL for the perturbation
173  Dual_ dummy;
174  CtrlInc_ mdx(dx);
175  H.multiply(mdx, dummy);
176 
177 // print log
178  Log::info() << std::endl << "TLM Linear Approximation Test: done." << std::endl;
179 }
180 
181 // -----------------------------------------------------------------------------
182 
183 template<typename MODEL, typename OBS>
184 void Minimizer<MODEL, OBS>::tlmPropagTest(const eckit::Configuration & config,
185  const CtrlInc_ & dx) {
186 /* TL propagation test:
187  calculate and store:
188  M'(dx)
189  where dx comes from minimization stopped at requested iteration
190 
191  TO DO:
192  write to file depending on requriements */
193 
194  if (config.has("online diagnostics")) {
195  const eckit::LocalConfiguration onlineDiag(config, "online diagnostics");
196  bool runTlmPropagTest = onlineDiag.getBool("tlm propagation test", false);
197 
198  if (runTlmPropagTest) {
199  // print log
200  Log::info() << "TLM Propagation Test - starting: " << outerIteration_
201  << std::endl << std::endl;
202 
203  // construct propagation matrix
204  const H_ H(J_);
205 
206  // run TL for the input perturbation: M'(dx)
207  Dual_ dummy;
208  CtrlInc_ mdx(dx);
209  H.multiply(mdx, dummy);
210 
211  // print log
212  Log::info() << std::endl << "TLM Propagation Test: done." << std::endl;
213  }
214  }
215 }
216 
217 // -----------------------------------------------------------------------------
218 
219 template<typename MODEL, typename OBS>
221 /* TL Taylor test:
222  ||M(x + p dx) - M(x)|| / ||M'(p dx)|| -> 1
223  where p is a scalar damping factor and dx is initialized
224  using randomization */
225 
226  Log::info() << std::endl
227  << "TLM Taylor Test: " << outerIteration_
228  << std::endl;
229 
230 // initialize incremets
231  CtrlInc_ dx(J_.jb());
232 
233 // randomize increments
234  J_.jb().randomize(dx);
235 
236 // run NL for unperturbed initial condition: M(x)
238  ControlVariable<MODEL, OBS> mxx(J_.jb().getBackground());
239  J_.runNL(mxx, pp);
240 
241 // run TL for the un-damped perturbation: M'(dx)
242  Dual_ dummy;
243  CtrlInc_ mdx(dx);
244  H.multiply(mdx, dummy);
245 
246 // loop over decreasing increments
247  for (unsigned int jj = 0; jj < 14; ++jj) {
248  // ||p M'(dx)||
249  CtrlInc_ pmdx(mdx);
250  pmdx *= 1./pow(10.0, jj);
251  double denom = sqrt(dot_product(pmdx, pmdx));
252 
253  // run perturbed NL: M(x+pdx)
254  ControlVariable<MODEL, OBS> mpertxx(J_.jb().getBackground());
255  CtrlInc_ pdx(dx);
256  pdx *= 1./pow(10.0, jj);
257  J_.addIncrement(mpertxx, pdx);
258  J_.runNL(mpertxx, pp);
259 
260  // ||M(x+pdx) - M(x)||
261  CtrlInc_ diff_nl(mdx, false);
262  diff_nl.state().diff(mpertxx.state(), mxx.state());
263  double nom = sqrt(dot_product(diff_nl, diff_nl));
264 
265  // print results
266  Log::info() << std::endl
267  << "TLM Taylor test: p = " << std::setw(8) << 1./pow(10.0, jj)
268  << ", ||M(x) - M(x+p dx)|| / ||p M'(dx)|| = "
269  << util::full_precision(1.+std::abs(1.-nom/denom))
270  << std::endl;
271  }
272  Log::info() << std::endl;
273 }
274 
275 // -----------------------------------------------------------------------------
276 
277 template<typename MODEL, typename OBS>
279  const H_ & H) {
280 /* Perform the adjoint test of the linear model
281  <M dx1, dx2> - <dx1, Mt dx2>)/<M dx1, dx2>
282  where dx1 and dx2 are increments obtained through randomization */
283 
284 // Model adjoint test
285  CtrlInc_ dx1(J_.jb());
286  CtrlInc_ dx2(J_.jb());
287 
288 // randomize increments
289  J_.jb().randomize(dx1);
290  J_.jb().randomize(dx2);
291 
292 /* zero out modvar and obsvar
293  since only model is being tested
294 */
295  dx1.obsVar().zero();
296  dx1.modVar().zero();
297  dx2.obsVar().zero();
298  dx2.modVar().zero();
299 
300 // run TL
301  Dual_ dummy;
302  CtrlInc_ mdx1(dx1);
303  H.multiply(mdx1, dummy, false);
304 
305 // run ADJ
306  dummy.zero();
307  CtrlInc_ mtdx2(dx2);
308  mtdx2.state().updateTime(mdx1.state().validTime() - dx1.state().validTime());
309  Ht.multiply(dummy, mtdx2, false);
310 
311 // calculate FWD < M dx1, dx2 >
312  double adj_tst_fwd = dot_product(mdx1, dx2);
313 
314 // calculate BWD < dx1, Mt dx2 >
315  double adj_tst_bwd = dot_product(dx1, mtdx2);
316 
317 // print results
318  Log::info() << "Model Adjoint Test: " << outerIteration_ << std::endl
319  << util::PrintAdjTest(adj_tst_fwd, adj_tst_bwd, "M")
320  << std::endl << std::endl;
321 }
322 
323 // -----------------------------------------------------------------------------
324 
325 template<typename MODEL, typename OBS>
327  const H_ & H) {
328 /* Perform the adjoint test of the linear observation operator
329  (<H dx, dy> - <dx, Ht dy>)/<H dx, dy>
330  where dx is an increment obtained through randomization and
331  dy is a DualVector obtained by first creating an increment
332  through randomization and then applying a linear observation
333  operator */
334 
335 // Model adjoint test
336  CtrlInc_ dx1(J_.jb());
337  CtrlInc_ dx2(J_.jb());
338  CtrlInc_ hthdx2(J_.jb());
339  Dual_ hdx1;
340  Dual_ hdx2;
341 
342 // randomize increments
343  J_.jb().randomize(dx1);
344  J_.jb().randomize(dx2);
345 
346 // run TL
347  H.multiply(dx1, hdx1, true);
348  H.multiply(dx2, hdx2, true);
349 
350 // run ADJ
351  J_.zeroAD(hthdx2);
352  Ht.multiply(hdx2, hthdx2, true);
353 
354 // calculate FWD < H dx1, hdx2 >
355  double adj_tst_fwd = dot_product(hdx1, hdx2);
356 
357 // calculate BWD < dx1, Ht hdx2 >
358  double adj_tst_bwd = dot_product(dx1, hthdx2);
359 
360 // print results
361  Log::info() << "Obs Adjoint Test: " << outerIteration_ << std::endl
362  << util::PrintAdjTest(adj_tst_fwd, adj_tst_bwd, "H")
363  << std::endl << std::endl;
364 }
365 
366 // -----------------------------------------------------------------------------
367 
368 /// Minimizer Factory
369 template <typename MODEL, typename OBS>
370 class MinFactory {
372  public:
373  static Minimizer<MODEL, OBS> * create(const eckit::Configuration &, const CostFct_ &);
374  virtual ~MinFactory() = default;
375  protected:
376  explicit MinFactory(const std::string &);
377  private:
378  virtual Minimizer<MODEL, OBS> * make(const eckit::Configuration &, const CostFct_ &) = 0;
379  static std::map < std::string, MinFactory<MODEL, OBS> * > & getMakers() {
380  static std::map < std::string, MinFactory<MODEL, OBS> * > makers_;
381  return makers_;
382  }
383 };
384 
385 // -----------------------------------------------------------------------------
386 
387 template<class MODEL, class OBS, class FCT>
388 class MinMaker : public MinFactory<MODEL, OBS> {
390  virtual Minimizer<MODEL, OBS> * make(const eckit::Configuration & conf, const CostFct_ & J) {
391  return new FCT(conf, J);
392  }
393  public:
394  explicit MinMaker(const std::string & name) : MinFactory<MODEL, OBS>(name) {}
395 };
396 
397 // =============================================================================
398 
399 template <typename MODEL, typename OBS>
400 MinFactory<MODEL, OBS>::MinFactory(const std::string & name) {
401  if (getMakers().find(name) != getMakers().end()) {
402  throw std::runtime_error(name + " already registered in minimizer factory.");
403  }
404  getMakers()[name] = this;
405 }
406 
407 // -----------------------------------------------------------------------------
408 
409 template <typename MODEL, typename OBS>
410 Minimizer<MODEL, OBS>* MinFactory<MODEL, OBS>::create(const eckit::Configuration & config,
411  const CostFct_ & J) {
412  std::string id = config.getString("algorithm");
413  Log::info() << "Minimizer algorithm=" << id << std::endl;
414  typename std::map<std::string, MinFactory<MODEL, OBS>*>::iterator j = getMakers().find(id);
415  if (j == getMakers().end()) {
416  throw std::runtime_error(id + " does not exist in minimizer factory.");
417  }
418  return (*j).second->make(config, J);
419 }
420 
421 // -----------------------------------------------------------------------------
422 
423 } // namespace oops
424 
425 #endif // OOPS_ASSIMILATION_MINIMIZER_H_
ObsAuxIncrs_ & obsVar()
Get augmented observation control variable.
Increment_ & state()
Get state control variable.
ModelAuxIncr_ & modVar()
Get augmented model control variable.
Control variable.
State_ & state()
Get state control variable.
Cost Function.
Definition: CostFunction.h:53
Container of dual space vectors for all terms of the cost function.
Definition: DualVector.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
Minimizer Factory.
Definition: Minimizer.h:370
virtual ~MinFactory()=default
static std::map< std::string, MinFactory< MODEL, OBS > * > & getMakers()
Definition: Minimizer.h:379
MinFactory(const std::string &)
Definition: Minimizer.h:400
static Minimizer< MODEL, OBS > * create(const eckit::Configuration &, const CostFct_ &)
Definition: Minimizer.h:410
virtual Minimizer< MODEL, OBS > * make(const eckit::Configuration &, const CostFct_ &)=0
CostFunction< MODEL, OBS > CostFct_
Definition: Minimizer.h:371
MinMaker(const std::string &name)
Definition: Minimizer.h:394
virtual Minimizer< MODEL, OBS > * make(const eckit::Configuration &conf, const CostFct_ &J)
Definition: Minimizer.h:390
CostFunction< MODEL, OBS > CostFct_
Definition: Minimizer.h:389
A Minimizer knows how to minimize a cost function.
Definition: Minimizer.h:37
HtMatrix< MODEL, OBS > Ht_
Definition: Minimizer.h:42
DualVector< MODEL, OBS > Dual_
Definition: Minimizer.h:40
Minimizer(const CostFct_ &J)
Definition: Minimizer.h:46
virtual const std::string classname() const =0
virtual ~Minimizer()
Definition: Minimizer.h:47
void tlmTests(const eckit::Configuration &)
Definition: Minimizer.h:97
State< MODEL > State_
Definition: Minimizer.h:43
ControlIncrement< MODEL, OBS > * minimize(const eckit::Configuration &)
Definition: Minimizer.h:72
void adjObsTest(const Ht_ &, const H_ &)
Definition: Minimizer.h:326
int outerIteration_
Definition: Minimizer.h:65
void tlmTaylorTest(const H_ &)
Definition: Minimizer.h:220
void tlmPropagTest(const eckit::Configuration &config, const CtrlInc_ &)
Definition: Minimizer.h:184
void adjModelTest(const Ht_ &, const H_ &)
Definition: Minimizer.h:278
void adjTests(const eckit::Configuration &)
Definition: Minimizer.h:120
const CostFct_ & J_
Definition: Minimizer.h:64
virtual ControlIncrement< MODEL, OBS > * doMinimize(const eckit::Configuration &)=0
void tlmApproxTest(const H_ &)
Definition: Minimizer.h:143
HMatrix< MODEL, OBS > H_
Definition: Minimizer.h:41
ControlIncrement< MODEL, OBS > CtrlInc_
Definition: Minimizer.h:39
CostFunction< MODEL, OBS > CostFct_
Definition: Minimizer.h:38
void zero()
Zero out this ModelAuxIncrement.
Control model post processing.
Definition: PostProcessor.h:30
State class used in oops; subclass of interface class interface::State.
void updateTime(const util::Duration &dt)
Updates this Increment's valid time by dt (used in PseudoModel)
void diff(const State_ &state1, const State_ &state2)
Set this Increment to be difference between state1 and state2.
const util::DateTime validTime() const
Accessor to the time of this Increment.
The namespace for the main oops code.
void writeIncrement(const eckit::Configuration &config, const ControlIncrement< MODEL, OBS > &dx, const int &loop)