11 #ifndef OOPS_ASSIMILATION_RPLANCZOSMINIMIZER_H_
12 #define OOPS_ASSIMILATION_RPLANCZOSMINIMIZER_H_
25 #include "oops/util/dot_product.h"
26 #include "oops/util/formats.h"
27 #include "oops/util/Logger.h"
76 const std::string
classname()
const override {
return "RPLanczosMinimizer";}
83 const int &,
const double &,
Dual_ &,
const double &)
override;
88 template<
typename MODEL,
typename OBS>
91 const int & maxiter,
const double & tolerance,
92 Dual_ & dy,
const double & sigma) {
109 std::vector<Dual_> vVEC;
110 std::vector<Dual_> tVEC;
111 std::vector<Dual_> zVEC;
113 vVEC.reserve(maxiter+1);
114 tVEC.reserve(maxiter+1);
115 zVEC.reserve(maxiter+1);
117 std::vector<double> vpVEC;
118 std::vector<double> tpVEC;
119 std::vector<double> zpVEC;
120 std::vector<double> alphas;
121 std::vector<double> betas;
122 std::vector<double> y;
123 std::vector<double> d;
132 ttp = dot_product(dy, zz) + sigma*zzp;
134 double normReduction = 1.0;
135 double beta0 = sqrt(dot_product(rr, tt) + rrp*ttp);
162 zpVEC.push_back(zzp);
163 tpVEC.push_back(ttp);
166 Log::info() << std::endl;
167 while (jiter < maxiter) {
168 Log::info() <<
"RPLanczos Starting Iteration " << jiter+1 << std::endl;
174 ww.
axpy(-beta, vold);
177 double alpha = dot_product(tt, ww) + ttp*wwp;
183 for (
int iiter = 0; iiter < jiter; ++iiter) {
184 double proj = dot_product(ww, tVEC[iiter]) + wwp * tpVEC[iiter];
185 ww.
axpy(-proj, vVEC[iiter]);
186 wwp -= proj * vpVEC[iiter];
196 ttp = dot_product(dy, zz) + sigma*zzp;
198 beta = sqrt(dot_product(tt, ww) + ttp*wwp);
215 zpVEC.push_back(zzp);
216 tpVEC.push_back(ttp);
218 alphas.push_back(alpha);
222 y.push_back(beta0/alpha);
226 for (
int iiter = 0; iiter <= jiter; ++iiter) {
227 d.push_back(beta0*(dot_product(tVEC[0], vVEC[iiter]) + tpVEC[0]*vpVEC[iiter]));
233 double rznorm = beta*std::abs(y[jiter]);
235 normReduction = rznorm/beta0;
237 betas.push_back(beta);
239 Log::info() <<
"RPLanczos end of iteration " << jiter+1 << std::endl;
244 if (normReduction < tolerance) {
245 Log::info() <<
"RPLanczos: Achieved required reduction in residual norm." << std::endl;
251 for (
int iiter = 0; iiter < jiter; ++iiter) {
252 vv.
axpy(y[iiter], zVEC[iiter]);
253 vvp += y[iiter]*zpVEC[iiter];
256 Log::info() <<
"RPLanczos: end" << std::endl;
258 return normReduction;
Container of dual space vectors for all terms of the cost function.
void axpy(const double, const DualVector &)
void multiply(const Dual_ &dy, Dual_ &dz) const
void multiply(const VECTOR &a, VECTOR &b) const
double solve(Dual_ &, double &, Dual_ &, const HBHt_ &, const Rinv_ &, const int &, const double &, Dual_ &, const double &) override
RPLanczosMinimizer(const eckit::Configuration &, const CostFct_ &J)
RinvMatrix< MODEL, OBS > Rinv_
CostFunction< MODEL, OBS > CostFct_
const std::string classname() const override
HBHtMatrix< MODEL, OBS > HBHt_
DualVector< MODEL, OBS > Dual_
void multiply(const Dual_ &dx, Dual_ &dy) const
The namespace for the main oops code.
void TriDiagSolve(const std::vector< double > &diag, const std::vector< double > &sub, const std::vector< double > &rhs, std::vector< double > &sol)
void printNormReduction(int iteration, const double &grad, const double &norm)