11 #ifndef OOPS_ASSIMILATION_MINIMIZER_H_
12 #define OOPS_ASSIMILATION_MINIMIZER_H_
16 #include <boost/noncopyable.hpp>
18 #include "eckit/config/Configuration.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"
37 template<
typename MODEL,
typename OBS>
class Minimizer :
private boost::noncopyable {
54 void adjTests(
const eckit::Configuration &);
58 void tlmTests(
const eckit::Configuration &);
70 template<
typename MODEL,
typename OBS>
74 this->tlmTests(config);
77 this->adjTests(config);
86 this->tlmPropagTest(config, *dx);
96 template<
typename MODEL,
typename OBS>
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);
110 if (outerIteration_ == 0 && runTlmApproxTest) this->tlmApproxTest(H);
113 if (outerIteration_ == 0 && runTlmTaylorTest) this->tlmTaylorTest(H);
119 template<
typename MODEL,
typename OBS>
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);
133 if (runAdjTlmTest) this->adjModelTest(Ht, H);
136 if (runAdjObsTest) this->adjObsTest(Ht, H);
142 template<
typename MODEL,
typename OBS>
152 Log::info() << std::endl
153 <<
"TLM Linear Approximation Test - starting: " << outerIteration_
154 << std::endl << std::endl;
160 J_.jb().randomize(dx);
169 J_.addIncrement(mpertxx, dx);
170 J_.runNL(mpertxx, pp);
178 Log::info() << std::endl <<
"TLM Linear Approximation Test: done." << std::endl;
183 template<
typename MODEL,
typename OBS>
194 if (config.has(
"online diagnostics")) {
195 const eckit::LocalConfiguration onlineDiag(config,
"online diagnostics");
196 bool runTlmPropagTest = onlineDiag.getBool(
"tlm propagation test",
false);
198 if (runTlmPropagTest) {
200 Log::info() <<
"TLM Propagation Test - starting: " << outerIteration_
201 << std::endl << std::endl;
212 Log::info() << std::endl <<
"TLM Propagation Test: done." << std::endl;
219 template<
typename MODEL,
typename OBS>
226 Log::info() << std::endl
227 <<
"TLM Taylor Test: " << outerIteration_
234 J_.jb().randomize(dx);
247 for (
unsigned int jj = 0; jj < 14; ++jj) {
250 pmdx *= 1./pow(10.0, jj);
251 double denom = sqrt(dot_product(pmdx, pmdx));
256 pdx *= 1./pow(10.0, jj);
257 J_.addIncrement(mpertxx, pdx);
258 J_.runNL(mpertxx, pp);
263 double nom = sqrt(dot_product(diff_nl, diff_nl));
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))
272 Log::info() << std::endl;
277 template<
typename MODEL,
typename OBS>
289 J_.jb().randomize(dx1);
290 J_.jb().randomize(dx2);
312 double adj_tst_fwd = dot_product(mdx1, dx2);
315 double adj_tst_bwd = dot_product(dx1, mtdx2);
318 Log::info() <<
"Model Adjoint Test: " << outerIteration_ << std::endl
319 << util::PrintAdjTest(adj_tst_fwd, adj_tst_bwd,
"M")
320 << std::endl << std::endl;
325 template<
typename MODEL,
typename OBS>
343 J_.jb().randomize(dx1);
344 J_.jb().randomize(dx2);
355 double adj_tst_fwd = dot_product(hdx1, hdx2);
358 double adj_tst_bwd = dot_product(dx1, hthdx2);
361 Log::info() <<
"Obs Adjoint Test: " << outerIteration_ << std::endl
362 << util::PrintAdjTest(adj_tst_fwd, adj_tst_bwd,
"H")
363 << std::endl << std::endl;
369 template <
typename MODEL,
typename OBS>
379 static std::map < std::string, MinFactory<MODEL, OBS> * > &
getMakers() {
380 static std::map < std::string, MinFactory<MODEL, OBS> * > makers_;
387 template<
class MODEL,
class OBS,
class FCT>
391 return new FCT(conf, J);
399 template <
typename MODEL,
typename OBS>
401 if (getMakers().find(name) != getMakers().end()) {
402 throw std::runtime_error(name +
" already registered in minimizer factory.");
404 getMakers()[name] =
this;
409 template <
typename MODEL,
typename OBS>
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.");
418 return (*j).second->make(config, J);
ObsAuxIncrs_ & obsVar()
Get augmented observation control variable.
Increment_ & state()
Get state control variable.
ModelAuxIncr_ & modVar()
Get augmented model control variable.
State_ & state()
Get state control variable.
Container of dual space vectors for all terms of the cost function.
void multiply(CtrlInc_ &dx, DualVector< MODEL, OBS > &dy, const bool idModel=false) const
void multiply(const DualVector< MODEL, OBS > &dy, ControlIncrement< MODEL, OBS > &dx, const bool idModel=false) const
virtual ~MinFactory()=default
static std::map< std::string, MinFactory< MODEL, OBS > * > & getMakers()
MinFactory(const std::string &)
static Minimizer< MODEL, OBS > * create(const eckit::Configuration &, const CostFct_ &)
virtual Minimizer< MODEL, OBS > * make(const eckit::Configuration &, const CostFct_ &)=0
CostFunction< MODEL, OBS > CostFct_
MinMaker(const std::string &name)
virtual Minimizer< MODEL, OBS > * make(const eckit::Configuration &conf, const CostFct_ &J)
CostFunction< MODEL, OBS > CostFct_
A Minimizer knows how to minimize a cost function.
HtMatrix< MODEL, OBS > Ht_
DualVector< MODEL, OBS > Dual_
Minimizer(const CostFct_ &J)
virtual const std::string classname() const =0
void tlmTests(const eckit::Configuration &)
ControlIncrement< MODEL, OBS > * minimize(const eckit::Configuration &)
void adjObsTest(const Ht_ &, const H_ &)
void tlmTaylorTest(const H_ &)
void tlmPropagTest(const eckit::Configuration &config, const CtrlInc_ &)
void adjModelTest(const Ht_ &, const H_ &)
void adjTests(const eckit::Configuration &)
virtual ControlIncrement< MODEL, OBS > * doMinimize(const eckit::Configuration &)=0
void tlmApproxTest(const H_ &)
ControlIncrement< MODEL, OBS > CtrlInc_
CostFunction< MODEL, OBS > CostFct_
void zero()
Zero out this ModelAuxIncrement.
Control model post processing.
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)