OOPS
FGMRESMinimizer.h
Go to the documentation of this file.
1 /*
2  * (C) Copyright 2009-2016 ECMWF.
3  *
4  * This software is licensed under the terms of the Apache Licence Version 2.0
5  * which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
6  * In applying this licence, ECMWF does not waive the privileges and immunities
7  * granted to it by virtue of its status as an intergovernmental organisation nor
8  * does it submit to any jurisdiction.
9  */
10 
11 #ifndef OOPS_ASSIMILATION_FGMRESMINIMIZER_H_
12 #define OOPS_ASSIMILATION_FGMRESMINIMIZER_H_
13 
14 #include <string>
15 
22 
23 namespace oops {
24 
25 /// FGMRES Minimizer
26 /*!
27  * Implements the Flexible GMRES solver for Ax=b.
28  */
29 
30 // -----------------------------------------------------------------------------
31 
32 template<typename MODEL, typename OBS> class FGMRESMinimizer : public PrimalMinimizer<MODEL, OBS> {
37 
38  public:
39  const std::string classname() const override {return "FGMRESMinimizer";}
40  FGMRESMinimizer(const eckit::Configuration &, const CostFct_ & J)
41  : PrimalMinimizer<MODEL, OBS>(J) {}
43 
44  private:
45  double solve(CtrlInc_ &, const CtrlInc_ &,
46  const Hessian_ &, const Bmat_ &,
47  const int, const double) override;
48 };
49 
50 // =============================================================================
51 
52 template<typename MODEL, typename OBS>
54  const Hessian_ & hessian, const Bmat_ & B,
55  const int ninner, const double gnreduc) {
56 // Solve the linear system
57  double reduc = FGMRES(dx, rhs, hessian, B, ninner, gnreduc);
58  return reduc;
59 }
60 
61 // -----------------------------------------------------------------------------
62 
63 } // namespace oops
64 
65 #endif // OOPS_ASSIMILATION_FGMRESMINIMIZER_H_
FGMRES solver for Ax=b.
The matrix.
Definition: BMatrix.h:27
Cost Function.
Definition: CostFunction.h:53
FGMRES Minimizer.
const std::string classname() const override
BMatrix< MODEL, OBS > Bmat_
double solve(CtrlInc_ &, const CtrlInc_ &, const Hessian_ &, const Bmat_ &, const int, const double) override
HessianMatrix< MODEL, OBS > Hessian_
FGMRESMinimizer(const eckit::Configuration &, const CostFct_ &J)
CostFunction< MODEL, OBS > CostFct_
ControlIncrement< MODEL, OBS > CtrlInc_
The Hessian matrix: .
Definition: HessianMatrix.h:34
Primal Minimizer.
The namespace for the main oops code.
double FGMRES(VECTOR &x, const VECTOR &b, const AMATRIX &A, const PMATRIX &precond, const int maxiter, const double tolerance)
Definition: FGMRES.h:69