FV3 Bundle
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"
27 #include "oops/interface/State.h"
28 #include "oops/util/abor1_cpp.h"
29 #include "oops/util/dot_product.h"
30 #include "oops/util/formats.h"
31 #include "oops/util/Logger.h"
32 #include "oops/util/PrintAdjTest.h"
33 
34 namespace oops {
35 
36 // -----------------------------------------------------------------------------
37 
38 /// A Minimizer knows how to minimize a cost function
39 template<typename MODEL> class Minimizer : private boost::noncopyable {
43  typedef HMatrix<MODEL> H_;
46 
47  public:
48  explicit Minimizer(const CostFct_ & J): J_(J), outerIteration_(0) {}
49  virtual ~Minimizer() {}
50  ControlIncrement<MODEL> * minimize(const eckit::Configuration &);
51  virtual const std::string classname() const = 0;
52 
53  private:
54  virtual ControlIncrement<MODEL> * doMinimize(const eckit::Configuration &) = 0;
55 
56  void adjTests(const eckit::Configuration &);
57  void adjModelTest(const Ht_ &, const H_ &);
58  void adjObsTest(const Ht_ &, const H_ &);
59 
60  void tlmTests(const eckit::Configuration &);
61  void tlmApproxTest(const H_ &);
62  void tlmTaylorTest(const H_ &);
63 
64  void tlmPropagTest(const eckit::Configuration & config, const CtrlInc_ &);
65 
66  const CostFct_ & J_;
68 };
69 
70 // -----------------------------------------------------------------------------
71 
72 template<typename MODEL>
73 ControlIncrement<MODEL> * Minimizer<MODEL>::minimize(const eckit::Configuration & config) {
74  // TLM tests
75  this->tlmTests(config);
76 
77  // ADJ tests
78  this->adjTests(config);
79 
80  // Minimize
81  ControlIncrement<MODEL> * dx = this->doMinimize(config);
82 
83  // TLM propagation test
84  this->tlmPropagTest(config, *dx);
85 
86  // Update outer loop counter
87  outerIteration_++;
88 
89  return dx;
90 }
91 
92 // -----------------------------------------------------------------------------
93 
94 template<typename MODEL>
95 void Minimizer<MODEL>::tlmTests(const eckit::Configuration & config) {
96 // Tangent Linear tests
97 
98  if (config.has("onlineDiagnostics")) {
99  const eckit::LocalConfiguration onlineDiag(config, "onlineDiagnostics");
100  bool runTlmTaylorTest = onlineDiag.getBool("tlmTaylorTest");
101  bool runTlmApproxTest = onlineDiag.getBool("tlmApproxTest");
102 
103  const H_ H(J_);
104  const Ht_ Ht(J_);
105 
106  // Online tests
107  // TLM linear approximation test
108  if (outerIteration_ == 0 && runTlmApproxTest) this->tlmApproxTest(H);
109 
110  // TLM Taylor test
111  if (outerIteration_ == 0 && runTlmTaylorTest) this->tlmTaylorTest(H);
112  }
113 }
114 
115 // -----------------------------------------------------------------------------
116 
117 template<typename MODEL>
118 void Minimizer<MODEL>::adjTests(const eckit::Configuration & config) {
119 // Adjoint tests
120 
121  if (config.has("onlineDiagnostics")) {
122  const eckit::LocalConfiguration onlineDiag(config, "onlineDiagnostics");
123  bool runAdjTlmTest = onlineDiag.getBool("adjTlmTest");
124  bool runAdjObsTest = onlineDiag.getBool("adjObsTest");
125 
126  const H_ H(J_);
127  const Ht_ Ht(J_);
128 
129  // Online tests
130  // Model adjoint test
131  if (runAdjTlmTest) this->adjModelTest(Ht, H);
132 
133  // Obs adjoint test
134  if (runAdjObsTest) this->adjObsTest(Ht, H);
135  }
136 }
137 
138 // -----------------------------------------------------------------------------
139 
140 template<typename MODEL>
142 /* TL approx test:
143  calculate and store:
144  M(x + dx), M(x), M'(dx)
145  where dx is initialized using randomization
146 
147  TO DO:
148  write to file depending on requriements */
149 
150  Log::info() << std::endl
151  << "TLM Linear Approximation Test - starting: " << outerIteration_
152  << std::endl << std::endl;
153 
154 // initialize incremets
155  CtrlInc_ dx(J_.jb());
156 
157 // randomize increments
158  J_.jb().randomize(dx);
159 
160 // run NL for unperturbed initial condition
162  ControlVariable<MODEL> mxx(J_.jb().getBackground());
163  J_.runNL(mxx, pp);
164 
165 // run NL for perturbed initial condition
166  ControlVariable<MODEL> mpertxx(J_.jb().getBackground());
167  J_.addIncrement(mpertxx, dx);
168  J_.runNL(mpertxx, pp);
169 
170 // run TL for the perturbation
171  Dual_ dummy;
172  CtrlInc_ mdx(dx);
173  H.multiply(mdx, dummy);
174 
175 // print log
176  Log::info() << std::endl << "TLM Linear Approximation Test: done." << std::endl;
177 }
178 
179 // -----------------------------------------------------------------------------
180 
181 template<typename MODEL>
182 void Minimizer<MODEL>::tlmPropagTest(const eckit::Configuration & config,
183  const CtrlInc_ & dx) {
184 /* TL propagation test:
185  calculate and store:
186  M'(dx)
187  where dx comes from minimization stopped at requested iteration
188 
189  TO DO:
190  write to file depending on requriements */
191 
192  if (config.has("onlineDiagnostics")) {
193  const eckit::LocalConfiguration onlineDiag(config, "onlineDiagnostics");
194  bool runTlmPropagTest = onlineDiag.getBool("tlmPropagTest");
195 
196  if (runTlmPropagTest) {
197  // print log
198  Log::info() << "TLM Propagation Test - starting: " << outerIteration_
199  << std::endl << std::endl;
200 
201  // construct propagation matrix
202  const H_ H(J_);
203 
204  // run TL for the input perturbation: M'(dx)
205  Dual_ dummy;
206  CtrlInc_ mdx(dx);
207  H.multiply(mdx, dummy);
208 
209  // print log
210  Log::info() << std::endl << "TLM Propagation Test: done." << std::endl;
211  }
212  }
213 }
214 
215 // -----------------------------------------------------------------------------
216 
217 template<typename MODEL>
219 /* TL Taylor test:
220  ||M(x + p dx) - M(x)|| / ||M'(p dx)|| -> 1
221  where p is a scalar damping factor and dx is initialized
222  using randomization */
223 
224  Log::info() << std::endl
225  << "TLM Taylor Test: " << outerIteration_
226  << std::endl;
227 
228 // initialize incremets
229  CtrlInc_ dx(J_.jb());
230 
231 // randomize increments
232  J_.jb().randomize(dx);
233 
234 // run NL for unperturbed initial condition: M(x)
236  ControlVariable<MODEL> mxx(J_.jb().getBackground());
237  J_.runNL(mxx, pp);
238 
239 // run TL for the un-damped perturbation: M'(dx)
240  Dual_ dummy;
241  CtrlInc_ mdx(dx);
242  H.multiply(mdx, dummy);
243 
244 // loop over decreasing increments
245  for (unsigned int jj = 0; jj < 14; ++jj) {
246  // ||p M'(dx)||
247  CtrlInc_ pmdx(mdx);
248  pmdx *= 1./pow(10.0, jj);
249  double denom = sqrt(dot_product(pmdx, pmdx));
250 
251  // run perturbed NL: M(x+pdx)
252  ControlVariable<MODEL> mpertxx(J_.jb().getBackground());
253  CtrlInc_ pdx(dx);
254  pdx *= 1./pow(10.0, jj);
255  J_.addIncrement(mpertxx, pdx);
256  J_.runNL(mpertxx, pp);
257 
258  // ||M(x+pdx) - M(x)||
259  CtrlInc_ diff_nl(mdx, false);
260  diff_nl.state()[0].diff(mpertxx.state()[0], mxx.state()[0]);
261  double nom = sqrt(dot_product(diff_nl, diff_nl));
262 
263  // print results
264  Log::info() << std::endl
265  << "TLM Taylor test: p = " << std::setw(8) << 1./pow(10.0, jj)
266  << ", ||M(x) - M(x+p dx)|| / ||p M'(dx)|| = "
267  << util::full_precision(1.+std::abs(1.-nom/denom))
268  << std::endl;
269  }
270  Log::info() << std::endl;
271 }
272 
273 // -----------------------------------------------------------------------------
274 
275 template<typename MODEL>
277  const H_ & H) {
278 /* Perform the adjoint test of the linear model
279  <M dx1, dx2> - <dx1, Mt dx2>)/<M dx1, dx2>
280  where dx1 and dx2 are increments obtained through randomization */
281 
282 // Model adjoint test
283  CtrlInc_ dx1(J_.jb());
284  CtrlInc_ dx2(J_.jb());
285 
286 // randomize increments
287  J_.jb().randomize(dx1);
288  J_.jb().randomize(dx2);
289 
290 // run TL
291  Dual_ dummy;
292  CtrlInc_ mdx1(dx1);
293  H.multiply(mdx1, dummy, false);
294 
295 // run ADJ
296  dummy.zero();
297  CtrlInc_ mtdx2(dx2);
298  mtdx2.state()[0].updateTime(mdx1.state()[0].validTime() - dx1.state()[0].validTime());
299  Ht.multiply(dummy, mtdx2, false);
300 
301 // calculate FWD < M dx1, dx2 >
302  double adj_tst_fwd = dot_product(mdx1, dx2);
303 
304 // calculate BWD < dx1, Mt dx2 >
305  double adj_tst_bwd = dot_product(dx1, mtdx2);
306 
307 // print results
308  Log::info() << "Model Adjoint Test: " << outerIteration_ << std::endl
309  << util::PrintAdjTest(adj_tst_fwd, adj_tst_bwd, "M")
310  << std::endl << std::endl;
311 }
312 
313 // -----------------------------------------------------------------------------
314 
315 template<typename MODEL>
317  const H_ & H) {
318 /* Perform the adjoint test of the linear observation operator
319  (<H dx, dy> - <dx, Ht dy>)/<H dx, dy>
320  where dx is an increment obtained through randomization and
321  dy is a DualVector obtained by first creating an increment
322  through randomization and then applying a linear observation
323  operator */
324 
325 // Model adjoint test
326  CtrlInc_ dx1(J_.jb());
327  CtrlInc_ dx2(J_.jb());
328  CtrlInc_ hthdx2(J_.jb());
329  Dual_ hdx1;
330  Dual_ hdx2;
331 
332 // randomize increments
333  J_.jb().randomize(dx1);
334  J_.jb().randomize(dx2);
335 
336 // run TL
337  H.multiply(dx1, hdx1, true);
338  H.multiply(dx2, hdx2, true);
339 
340 // run ADJ
341  J_.zeroAD(hthdx2);
342  Ht.multiply(hdx2, hthdx2, true);
343 
344 // calculate FWD < H dx1, hdx2 >
345  double adj_tst_fwd = dot_product(hdx1, hdx2);
346 
347 // calculate BWD < dx1, Ht hdx2 >
348  double adj_tst_bwd = dot_product(dx1, hthdx2);
349 
350 // print results
351  Log::info() << "Obs Adjoint Test: " << outerIteration_ << std::endl
352  << util::PrintAdjTest(adj_tst_fwd, adj_tst_bwd, "H")
353  << std::endl << std::endl;
354 }
355 
356 // -----------------------------------------------------------------------------
357 
358 /// Minimizer Factory
359 template <typename MODEL>
360 class MinFactory {
362  public:
363  static Minimizer<MODEL> * create(const eckit::Configuration &, const CostFct_ &);
364  virtual ~MinFactory() { getMakers().clear(); }
365  protected:
366  explicit MinFactory(const std::string &);
367  private:
368  virtual Minimizer<MODEL> * make(const eckit::Configuration &, const CostFct_ &) = 0;
369  static std::map < std::string, MinFactory<MODEL> * > & getMakers() {
370  static std::map < std::string, MinFactory<MODEL> * > makers_;
371  return makers_;
372  }
373 };
374 
375 // -----------------------------------------------------------------------------
376 
377 template<class MODEL, class FCT>
378 class MinMaker : public MinFactory<MODEL> {
380  virtual Minimizer<MODEL> * make(const eckit::Configuration & conf, const CostFct_ & J) {
381  return new FCT(conf, J);
382  }
383  public:
384  explicit MinMaker(const std::string & name) : MinFactory<MODEL>(name) {}
385 };
386 
387 // =============================================================================
388 
389 template <typename MODEL>
390 MinFactory<MODEL>::MinFactory(const std::string & name) {
391  if (getMakers().find(name) != getMakers().end()) {
392  Log::error() << name << " already registered in minimizer factory." << std::endl;
393  ABORT("Element already registered in MinFactory.");
394  }
395  getMakers()[name] = this;
396 }
397 
398 // -----------------------------------------------------------------------------
399 
400 template <typename MODEL>
401 Minimizer<MODEL>* MinFactory<MODEL>::create(const eckit::Configuration & config,
402  const CostFct_ & J) {
403  std::string id = config.getString("algorithm");
404  Log::info() << "Minimizer algorithm=" << id << std::endl;
405  typename std::map<std::string, MinFactory<MODEL>*>::iterator j = getMakers().find(id);
406  if (j == getMakers().end()) {
407  Log::error() << id << " does not exist in minimizer factory." << std::endl;
408  ABORT("Element does not exist in MinFactory.");
409  }
410  return (*j).second->make(config, J);
411 }
412 
413 // -----------------------------------------------------------------------------
414 
415 } // namespace oops
416 
417 #endif // OOPS_ASSIMILATION_MINIMIZER_H_
The matrix.
Definition: HMatrix.h:33
Container of dual space vectors for all terms of the cost function.
Definition: DualVector.h:34
Minimizer Factory.
Definition: Minimizer.h:360
void tlmApproxTest(const H_ &)
Definition: Minimizer.h:141
ControlIncrement< MODEL > * minimize(const eckit::Configuration &)
Definition: Minimizer.h:73
CostFunction< MODEL > CostFct_
Definition: Minimizer.h:40
virtual Minimizer< MODEL > * make(const eckit::Configuration &, const CostFct_ &)=0
virtual ~Minimizer()
Definition: Minimizer.h:49
MinMaker(const std::string &name)
Definition: Minimizer.h:384
************************************************************************GNU Lesser General Public License **This file is part of the GFDL Flexible Modeling System(FMS). ! *! *FMS is free software without even the implied warranty of MERCHANTABILITY or *FITNESS FOR A PARTICULAR PURPOSE See the GNU General Public License *for more details **You should have received a copy of the GNU Lesser General Public *License along with FMS If see< http:! ***********************************************************************!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! !! MPP_TRANSMIT !! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! subroutine MPP_TRANSMIT_(put_data, put_len, to_pe, get_data, get_len, from_pe, block, tag, recv_request, send_request)!a message-passing routine intended to be reminiscent equally of both MPI and SHMEM!put_data and get_data are contiguous MPP_TYPE_ arrays!at each call, your put_data array is put to to_pe 's get_data! your get_data array is got from from_pe 's put_data!i.e we assume that typically(e.g updating halo regions) each PE performs a put _and_ a get!special PE designations:! NULL_PE:to disable a put or a get(e.g at boundaries)! ANY_PE:if remote PE for the put or get is to be unspecific! ALL_PES:broadcast and collect operations(collect not yet implemented)!ideally we would not pass length, but this f77-style call performs better(arrays passed by address, not descriptor)!further, this permits< length > contiguous words from an array of any rank to be passed(avoiding f90 rank conformance check)!caller is responsible for completion checks(mpp_sync_self) before and after integer, intent(in) ::put_len, to_pe, get_len, from_pe MPP_TYPE_, intent(in) ::put_data(*) MPP_TYPE_, intent(out) ::get_data(*) logical, intent(in), optional ::block integer, intent(in), optional ::tag integer, intent(out), optional ::recv_request, send_request logical ::block_comm integer ::i MPP_TYPE_, allocatable, save ::local_data(:) !local copy used by non-parallel code(no SHMEM or MPI) integer ::comm_tag integer ::rsize if(.NOT.module_is_initialized) call mpp_error(FATAL, 'MPP_TRANSMIT:You must first call mpp_init.') if(to_pe.EQ.NULL_PE .AND. from_pe.EQ.NULL_PE) return block_comm=.true. if(PRESENT(block)) block_comm=block if(debug) then call SYSTEM_CLOCK(tick) write(stdout_unit,'(a, i18, a, i6, a, 2i6, 2i8)')&'T=', tick, ' PE=', pe, ' MPP_TRANSMIT begin:to_pe, from_pe, put_len, get_len=', to_pe, from_pe, put_len, get_len end if comm_tag=DEFAULT_TAG if(present(tag)) comm_tag=tag!do put first and then get if(to_pe.GE.0 .AND. to_pe.LT.npes) then!use non-blocking sends if(debug .and.(current_clock.NE.0)) call SYSTEM_CLOCK(start_tick)!z1l:truly non-blocking send.! if(request(to_pe).NE.MPI_REQUEST_NULL) then !only one message from pe-> to_pe in queue *PE waiting for to_pe ! call error else get_len so only do gets but you cannot have a pure get with MPI call a get means do a wait to ensure put on remote PE is complete error call increase mpp_nml request_multiply call MPP_TRANSMIT end
Definition: conf.py:1
CostFunction< MODEL > CostFct_
Definition: Minimizer.h:379
HtMatrix< MODEL > Ht_
Definition: Minimizer.h:44
void diff(const State4D_ &, const State4D_ &)
Linear algebra operators.
Definition: Increment4D.h:171
void adjModelTest(const Ht_ &, const H_ &)
Definition: Minimizer.h:276
Cost Function.
Definition: CostFunction.h:56
virtual ~MinFactory()
Definition: Minimizer.h:364
Encapsulates the model state.
void multiply(CtrlInc_ &dx, DualVector< MODEL > &dy, const bool idModel=false) const
Definition: HMatrix.h:41
l_size ! loop over number of fields ke do j
Control variable.
void adjTests(const eckit::Configuration &)
Definition: Minimizer.h:118
character(len=32) name
static Minimizer< MODEL > * create(const eckit::Configuration &, const CostFct_ &)
Definition: Minimizer.h:401
virtual ControlIncrement< MODEL > * doMinimize(const eckit::Configuration &)=0
The namespace for the main oops code.
virtual const std::string classname() const =0
void adjObsTest(const Ht_ &, const H_ &)
Definition: Minimizer.h:316
Increment4D_ & state()
Get state control variable.
void multiply(const DualVector< MODEL > &dy, ControlIncrement< MODEL > &dx, const bool idModel=false) const
Definition: HtMatrix.h:40
Minimizer(const CostFct_ &J)
Definition: Minimizer.h:48
subroutine, public info(self)
integer error
Definition: mpp.F90:1310
CostFunction< MODEL > CostFct_
Definition: Minimizer.h:361
void tlmTests(const eckit::Configuration &)
Definition: Minimizer.h:95
virtual Minimizer< MODEL > * make(const eckit::Configuration &conf, const CostFct_ &J)
Definition: Minimizer.h:380
DualVector< MODEL > Dual_
Definition: Minimizer.h:42
int outerIteration_
Definition: Minimizer.h:67
State< MODEL > State_
Definition: Minimizer.h:45
static std::map< std::string, MinFactory< MODEL > *> & getMakers()
Definition: Minimizer.h:369
HMatrix< MODEL > H_
Definition: Minimizer.h:43
void tlmTaylorTest(const H_ &)
Definition: Minimizer.h:218
ControlIncrement< MODEL > CtrlInc_
Definition: Minimizer.h:41
Control model post processing.
Definition: PostProcessor.h:31
const CostFct_ & J_
Definition: Minimizer.h:66
A Minimizer knows how to minimize a cost function.
Definition: Minimizer.h:39
void tlmPropagTest(const eckit::Configuration &config, const CtrlInc_ &)
Definition: Minimizer.h:182
MinFactory(const std::string &)
Definition: Minimizer.h:390
The matrix.
Definition: HtMatrix.h:33