11 #ifndef OOPS_ASSIMILATION_SADDLEPOINTMINIMIZER_H_
12 #define OOPS_ASSIMILATION_SADDLEPOINTMINIMIZER_H_
19 #include "eckit/config/Configuration.h"
32 #include "oops/util/Logger.h"
51 const std::string
classname()
const override {
return "SaddlePointMinimizer";}
62 std::vector< SaddlePointVector<MODEL, OBS> >
xyVEC_;
63 std::vector< SaddlePointVector<MODEL, OBS> >
pqVEC_;
68 template<
typename MODEL,
typename OBS>
71 int ninner = config.getInt(
"ninner");
72 double gnreduc = config.getDouble(
"gradient norm reduction");
76 Log::info() <<
"SaddlePointMinimizer: max iter = " << ninner
77 <<
", requested norm reduction = " << gnreduc << std::endl;
82 for (
unsigned jj = 0; jj < J_.nterms(); ++jj) {
83 pdxx->
append(J_.jterm(jj).newDualVector());
92 for (
unsigned jj = 0; jj < J_.nterms(); ++jj) {
93 std::unique_ptr<GeneralizedDepartures> ww(J_.jterm(jj).newGradientFG());
94 pdfg->
append(J_.jterm(jj).multiplyCovar(*ww));
116 double reduc =
GMRESR(spdx, rhs, A, Pinv, ninner, gnreduc);
120 std::streamsize ss = Log::test().precision();
121 Log::test() <<
"SaddlePointMinimizer: reduction in residual norm = "
122 << std::setprecision(4) << reduc << std::setprecision(ss) << std::endl;
123 Log::info() <<
"SaddlePointMinimizer output" << *dx << std::endl;
Container of dual space vectors for all terms of the cost function.
void append(std::unique_ptr< GeneralizedDepartures > &&)
A Minimizer knows how to minimize a cost function.
DualVector< MODEL, OBS > Multipliers_
ControlIncrement< MODEL, OBS > CtrlInc_
CostFunction< MODEL, OBS > CostFct_
std::vector< SaddlePointVector< MODEL, OBS > > pqVEC_
SaddlePointMinimizer(const eckit::Configuration &, const CostFct_ &J)
std::vector< SaddlePointVector< MODEL, OBS > > xyVEC_
Minimizer< MODEL, OBS > Minimizer_
CtrlInc_ * doMinimize(const eckit::Configuration &) override
std::unique_ptr< CtrlInc_ > gradJb_
const std::string classname() const override
The preconditioner for the saddle-point minimizer.
Control vector for the saddle point formulation.
const CtrlInc_ & dx() const
Accessor method to get the dx_ component.
The namespace for the main oops code.
double GMRESR(VECTOR &xx, const VECTOR &bb, const AMATRIX &A, const PMATRIX &precond, const int maxiter, const double tolerance)
FullGMRES solver for Ax=b.