11 #ifndef OOPS_ASSIMILATION_DRPCGMINIMIZER_H_
12 #define OOPS_ASSIMILATION_DRPCGMINIMIZER_H_
18 #include "eckit/config/Configuration.h"
26 #include "oops/util/dot_product.h"
27 #include "oops/util/formats.h"
28 #include "oops/util/Logger.h"
86 const std::string
classname()
const override {
return "DRPCGMinimizer";}
92 const double,
const double,
const int,
const double)
override;
99 template<
typename MODEL,
typename OBS>
106 template<
typename MODEL,
typename OBS>
109 const double costJ0Jb,
const double costJ0JoJc,
110 const int maxiter,
const double tolerance) {
124 std::vector<CtrlInc_> vvecs;
125 std::vector<CtrlInc_> zvecs;
126 std::vector<double> scals;
128 vvecs.reserve(maxiter+1);
129 zvecs.reserve(maxiter+1);
130 scals.reserve(maxiter+1);
133 const double costJ0 = costJ0Jb + costJ0JoJc;
139 Log::info() <<
"normr0 " << dot_product(rr, rr) << std::endl;
142 lmp_.multiply(rr, pr);
150 double dotRr0 = dot_product(rr, zz);
151 double normReduction = 1.0;
152 double rdots = dotRr0;
153 double rdots_old = 0.0;
157 scals.push_back(1.0/dotRr0);
159 Log::info() << std::endl;
160 for (
int jiter = 0; jiter < maxiter; ++jiter) {
161 Log::info() <<
" DRPCG Starting Iteration " << jiter+1 << std::endl;
165 double beta = rdots/rdots_old;
181 double rho = dot_product(pp, qq);
182 double alpha = rdots/rho;
193 double costJ = costJ0 - 0.5 * dot_product(dx, r0);
195 double costJb = costJ0Jb + 0.5 * dot_product(dx, dxh);
197 double costJoJc = costJ - costJb;
200 for (
int jj = 0; jj < jiter; ++jj) {
201 double proj = scals[jj] * dot_product(rr, zvecs[jj]);
202 rr.
axpy(-proj, vvecs[jj]);
206 lmp_.multiply(rr, pr);
212 rdots = dot_product(rr, zz);
214 Log::info() <<
"rdots " << rdots <<
" iteration " << jiter << std::endl;
217 normReduction = sqrt(rdots/dotRr0);
219 Log::info() <<
"DRPCG end of iteration " << jiter+1 << std::endl;
224 lmp_.push(pp, hh, qq, rho);
226 if (normReduction < tolerance) {
227 Log::info() <<
"DRPCG: Achieved required reduction in residual norm." << std::endl;
233 scals.push_back(1.0/rdots);
239 return normReduction;
void multiply(const CtrlInc_ &x, CtrlInc_ &bx) const
void axpy(const double, const ControlIncrement &)
DR (Derber and Rosati) Minimizers.
QNewtonLMP< CtrlInc_, Bmat_ > lmp_
CostFunction< MODEL, OBS > CostFct_
ControlIncrement< MODEL, OBS > CtrlInc_
BMatrix< MODEL, OBS > Bmat_
HtRinvHMatrix< MODEL, OBS > HtRinvH_
const std::string classname() const override
DRPCGMinimizer(const eckit::Configuration &, const CostFct_ &)
double solve(CtrlInc_ &, CtrlInc_ &, CtrlInc_ &, const Bmat_ &, const HtRinvH_ &, const double, const double, const int, const double) override
void multiply(const CtrlInc_ &dx, CtrlInc_ &dz) const
The namespace for the main oops code.
void printQuadraticCostFunction(int iteration, const double &costJ, const double &costJb, const double &costJoJc)
void printNormReduction(int iteration, const double &grad, const double &norm)