11 #ifndef OOPS_ASSIMILATION_DRPLANCZOSMINIMIZER_H_
12 #define OOPS_ASSIMILATION_DRPLANCZOSMINIMIZER_H_
14 #include <Eigen/Dense>
29 #include "oops/util/dot_product.h"
30 #include "oops/util/formats.h"
31 #include "oops/util/Logger.h"
32 #include "oops/util/printRunStats.h"
83 const std::string
classname()
const override {
return "DRPLanczosMinimizer";}
89 const double,
const double,
const int,
const double)
override;
93 std::vector<std::unique_ptr<CtrlInc_>>
hvecs_;
94 std::vector<std::unique_ptr<CtrlInc_>>
vvecs_;
95 std::vector<std::unique_ptr<CtrlInc_>>
zvecs_;
105 template<
typename MODEL,
typename OBS>
108 :
DRMinimizer<MODEL, OBS>(J), lmp_(conf), hvecs_(), vvecs_(), zvecs_(), alphas_(),
109 betas_(), diagConf_(conf) {}
113 template<
typename MODEL,
typename OBS>
116 const double costJ0Jb,
const double costJ0JoJc,
117 const int maxiter,
const double tolerance) {
118 util::printRunStats(
"DRPLanczos start");
127 std::vector<double> ss;
128 std::vector<double> dd;
131 const double costJ0 = costJ0Jb + costJ0JoJc;
133 lmp_.update(vvecs_, hvecs_, zvecs_, alphas_, betas_);
136 lmp_.multiply(vv, pr);
140 double beta = sqrt(dot_product(zz, vv));
141 const double beta0 = beta;
151 hvecs_.emplace_back(std::unique_ptr<CtrlInc_>(
new CtrlInc_(pr)));
153 zvecs_.emplace_back(std::unique_ptr<CtrlInc_>(
new CtrlInc_(zz)));
155 vvecs_.emplace_back(std::unique_ptr<CtrlInc_>(
new CtrlInc_(vv)));
157 double normReduction = 1.0;
159 Log::info() << std::endl;
160 for (
int jiter = 0; jiter < maxiter; ++jiter) {
161 Log::info() <<
"DRPLanczos Starting Iteration " << jiter+1 << std::endl;
162 util::printRunStats(
"DRPLanczos iteration " + std::to_string(jiter+1));
169 vv.
axpy(-beta, *vvecs_[jiter-1]);
173 double alpha = dot_product(zz, vv);
176 vv.
axpy(-alpha, *vvecs_[jiter]);
179 for (
int jj = 0; jj < jiter; ++jj) {
180 double proj = dot_product(vv, *zvecs_[jj]);
181 vv.
axpy(-proj, *vvecs_[jj]);
185 lmp_.multiply(vv, pr);
189 beta = sqrt(dot_product(zz, vv));
199 hvecs_.emplace_back(std::unique_ptr<CtrlInc_>(
new CtrlInc_(pr)));
201 zvecs_.emplace_back(std::unique_ptr<CtrlInc_>(
new CtrlInc_(zz)));
203 vvecs_.emplace_back(std::unique_ptr<CtrlInc_>(
new CtrlInc_(vv)));
205 alphas_.push_back(alpha);
208 ss.push_back(beta0/alpha);
212 dd.push_back(beta0*dot_product(*zvecs_[0], vv));
216 betas_.push_back(beta);
221 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() <<
"DRPLanczos end of iteration " << jiter+1 << std::endl;
238 if (normReduction < tolerance) {
239 Log::info() <<
"DRPLanczos: Achieved required reduction in residual norm." << std::endl;
248 for (
unsigned int jj = 0; jj < ss.size(); ++jj) {
249 dx.
axpy(ss[jj], *zvecs_[jj]);
250 dxh.
axpy(ss[jj], *hvecs_[jj]);
254 writeEigenvectors(diagConf_, alphas_, betas_, dd, zvecs_, hvecs_, HtRinvH, pr, vv, zz);
256 util::printRunStats(
"DRPLanczos end");
257 return normReduction;
void multiply(const CtrlInc_ &x, CtrlInc_ &bx) const
void zero()
Linear algebra operators.
void axpy(const double, const ControlIncrement &)
DR (Derber and Rosati) Minimizers.
std::vector< std::unique_ptr< CtrlInc_ > > hvecs_
std::vector< double > alphas_
std::vector< double > betas_
double solve(CtrlInc_ &, CtrlInc_ &, CtrlInc_ &, const Bmat_ &, const HtRinvH_ &, const double, const double, const int, const double) override
BMatrix< MODEL, OBS > Bmat_
std::vector< std::unique_ptr< CtrlInc_ > > zvecs_
std::vector< std::unique_ptr< CtrlInc_ > > vvecs_
CostFunction< MODEL, OBS > CostFct_
const std::string classname() const override
DRPLanczosMinimizer(const eckit::Configuration &, const CostFct_ &)
HtRinvHMatrix< MODEL, OBS > HtRinvH_
eckit::LocalConfiguration diagConf_
SpectralLMP< CtrlInc_ > lmp_
ControlIncrement< MODEL, OBS > CtrlInc_
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 TriDiagSolve(const std::vector< double > &diag, const std::vector< double > &sub, const std::vector< double > &rhs, std::vector< double > &sol)
void writeEigenvectors(const eckit::Configuration &diagConf, const std::vector< double > &diag, const std::vector< double > &sub, const std::vector< double > &rhs, std::vector< std::unique_ptr< ControlIncrement< MODEL, OBS >>> &zvecs, std::vector< std::unique_ptr< ControlIncrement< MODEL, OBS >>> &hvecs, const HtRinvHMatrix< MODEL, OBS > &HtRinvH, ControlIncrement< MODEL, OBS > &temp, ControlIncrement< MODEL, OBS > &eigenv, ControlIncrement< MODEL, OBS > &eigenz)
void printNormReduction(int iteration, const double &grad, const double &norm)