11 #ifndef OOPS_ASSIMILATION_DRPFOMMINIMIZER_H_
12 #define OOPS_ASSIMILATION_DRPFOMMINIMIZER_H_
27 #include "oops/util/dot_product.h"
28 #include "oops/util/formats.h"
29 #include "oops/util/Logger.h"
74 const std::string
classname()
const override {
return "DRPFOMMinimizer";}
80 const double,
const double,
const int,
const double)
override;
85 std::vector<std::unique_ptr<CtrlInc_>>
hvecs_;
86 std::vector<std::unique_ptr<CtrlInc_>>
vvecs_;
87 std::vector<std::unique_ptr<CtrlInc_>>
zvecs_;
94 template<
typename MODEL,
typename OBS>
98 hvecs_(), vvecs_(), zvecs_(), alphas_(), betas_() {}
102 template<
typename MODEL,
typename OBS>
105 const double costJ0Jb,
const double costJ0JoJc,
106 const int maxiter,
const double tolerance) {
115 std::vector<double> ss;
116 std::vector<double> dd;
117 std::vector< std::vector<double> > Hess;
118 std::vector< std::vector<double> > UpHess;
121 const double costJ0 = costJ0Jb + costJ0JoJc;
135 double beta = sqrt(dot_product(zz, vv));
136 const double beta0 = beta;
146 hvecs_.emplace_back(std::unique_ptr<CtrlInc_>(
new CtrlInc_(pr)));
148 zvecs_.emplace_back(std::unique_ptr<CtrlInc_>(
new CtrlInc_(zz)));
150 vvecs_.emplace_back(std::unique_ptr<CtrlInc_>(
new CtrlInc_(vv)));
153 Hess.resize(maxiter);
154 for (
int ii = 0; ii <= maxiter-1; ii++) {
155 Hess[ii].resize(maxiter + 1);
156 for (
int jj = 0; jj <= maxiter; jj++) {
161 double normReduction = 1.0;
163 Log::info() << std::endl;
164 for (
int jiter = 0; jiter < maxiter; ++jiter) {
165 Log::info() <<
"DRPFOM Starting Iteration " << jiter+1 << std::endl;
171 for (
int jj = 0; jj <= jiter; ++jj) {
172 Hess[jiter][jj] = dot_product(*zvecs_[jj], vv);
173 vv.
axpy(-Hess[jiter][jj], *vvecs_[jj]);
182 beta = sqrt(dot_product(zz, vv));
184 Hess[jiter][jiter+1] = beta;
194 hvecs_.emplace_back(std::unique_ptr<CtrlInc_>(
new CtrlInc_(pr)));
196 zvecs_.emplace_back(std::unique_ptr<CtrlInc_>(
new CtrlInc_(zz)));
198 vvecs_.emplace_back(std::unique_ptr<CtrlInc_>(
new CtrlInc_(vv)));
201 ss.push_back(beta0/Hess[0][0]);
205 dd.push_back(beta0*dot_product(*zvecs_[0], vv));
207 UpHess.resize(jiter+1);
208 for (
int ii = 0; ii <= jiter; ii++) {
209 UpHess[ii].resize(jiter+1);
210 for (
int jj = 0; jj <= jiter; jj++) {
211 UpHess[ii][jj] = Hess[ii][jj];
217 betas_.push_back(beta);
222 double costJ = costJ0;
223 double costJb = costJ0Jb;
224 for (
int jj = 0; jj < jiter+1; ++jj) {
225 costJ -= 0.5 * ss[jj] * dot_product(*zvecs_[jj], rr);
226 costJb += 0.5 * ss[jj] * dot_product(*vvecs_[jj], *zvecs_[jj]) * ss[jj];
228 double costJoJc = costJ - costJb;
231 double rznorm = beta*std::abs(ss[jiter]);
232 normReduction = rznorm/beta0;
234 Log::info() <<
"DRPFOM end of iteration " << jiter+1 << std::endl;
238 if (normReduction < tolerance) {
239 Log::info() <<
"DRPFOM: Achieved required reduction in residual norm." << std::endl;
245 for (
unsigned int jj = 0; jj < ss.size(); ++jj) {
246 dx.
axpy(ss[jj], *zvecs_[jj]);
247 dxh.
axpy(ss[jj], *hvecs_[jj]);
250 return normReduction;
void multiply(const CtrlInc_ &x, CtrlInc_ &bx) const
void axpy(const double, const ControlIncrement &)
DR (Derber and Rosati) Minimizers.
SpectralLMP< CtrlInc_ > lmp_
HtRinvHMatrix< MODEL, OBS > HtRinvH_
std::vector< double > alphas_
double solve(CtrlInc_ &, CtrlInc_ &, CtrlInc_ &, const Bmat_ &, const HtRinvH_ &, const double, const double, const int, const double) override
CostFunction< MODEL, OBS > CostFct_
std::vector< std::unique_ptr< CtrlInc_ > > vvecs_
std::vector< std::unique_ptr< CtrlInc_ > > zvecs_
const std::string classname() const override
std::vector< std::unique_ptr< CtrlInc_ > > hvecs_
ControlIncrement< MODEL, OBS > CtrlInc_
DRPFOMMinimizer(const eckit::Configuration &, const CostFct_ &)
std::vector< double > betas_
BMatrix< MODEL, OBS > Bmat_
void multiply(const CtrlInc_ &dx, CtrlInc_ &dz) const
The solvers represent matrices as objects that implement a "multiply" method.
The namespace for the main oops code.
void printQuadraticCostFunction(int iteration, const double &costJ, const double &costJb, const double &costJoJc)
void UpHessSolve(std::vector< std::vector< double > > &UpHess, const std::vector< double > &rhs, std::vector< double > &sol)
void printNormReduction(int iteration, const double &grad, const double &norm)