11 #ifndef OOPS_ASSIMILATION_DUALMINIMIZER_H_
12 #define OOPS_ASSIMILATION_DUALMINIMIZER_H_
18 #include "eckit/config/Configuration.h"
28 #include "oops/util/dot_product.h"
29 #include "oops/util/Logger.h"
57 const int &,
const double &,
Dual_ &,
const double &) = 0;
65 template<
typename MODEL,
typename OBS>
68 int ninner = config.getInt(
"ninner");
69 double gnreduc = config.getDouble(
"gradient norm reduction");
71 bool runOnlineAdjTest = config.getBool(
"online diagnostics.online adj test",
false);
74 gradJb_.reset(
new CtrlInc_(J_.jb().resolution(), *gradJb_));
76 gradJb_.reset(
new CtrlInc_(J_.jb()));
79 Log::info() << std::endl;
80 Log::info() << classname() <<
": max iter = " << ninner
81 <<
", requested norm reduction = " << gnreduc << std::endl;
85 const HBHt_ HBHt(J_, runOnlineAdjTest);
94 for (
unsigned jj = 0; jj < J_.nterms(); ++jj) {
95 vv.
append(J_.jterm(jj).newDualVector());
101 for (
unsigned jj = 0; jj < J_.nterms(); ++jj) {
102 rr.
append(J_.jterm(jj).newGradientFG());
107 const double rnorm = dot_product(rr, rr);
108 const double epsilon = config.getDouble(
"epsilon", std::numeric_limits<double>::epsilon());
109 Log::info() <<
"Initial RHS squared norm = " << rnorm << std::endl;
110 if (rnorm < epsilon) {
111 Log::info() <<
"RHS smaller than " << epsilon <<
", returning." << std::endl;
120 CtrlInc_ hdxfg(J_.jb().getFirstGuess());
126 J_.jb().addGradientFG(g0, *gradJb_);
128 double sigma = dot_product(J_.jb().getFirstGuess(), g0);
131 double reduc = this->solve(vv, vvp, rr, HBHt, Rinv, ninner, gnreduc, dy, sigma);
133 Log::test() << classname() <<
": reduction in residual norm = " << reduc << std::endl;
143 Log::info() << classname() <<
": Estimated Final Jb = "
144 << 0.5 * dot_product(*dx, dh) << std::endl;
145 Log::info() << classname() <<
" output" << *dx << std::endl;
void multiply(const CtrlInc_ &x, CtrlInc_ &bx) const
void axpy(const double, const ControlIncrement &)
DualVector< MODEL, OBS > Dual_
const std::string classname() const override=0
CtrlInc_ * doMinimize(const eckit::Configuration &) override
RinvMatrix< MODEL, OBS > Rinv_
DualMinimizer(const CostFct_ &J)
std::unique_ptr< CtrlInc_ > gradJb_
Minimizer< MODEL, OBS > Minimizer_
CostFunction< MODEL, OBS > CostFct_
HBHtMatrix< MODEL, OBS > HBHt_
ControlIncrement< MODEL, OBS > CtrlInc_
BMatrix< MODEL, OBS > Bmat_
virtual double solve(Dual_ &, double &, Dual_ &, const HBHt_ &, const Rinv_ &, const int &, const double &, Dual_ &, const double &)=0
Container of dual space vectors for all terms of the cost function.
void append(std::unique_ptr< GeneralizedDepartures > &&)
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
A Minimizer knows how to minimize a cost function.
The namespace for the main oops code.