FV3 Bundle
GMRESRMinimizer.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_GMRESRMINIMIZER_H_
12 #define OOPS_ASSIMILATION_GMRESRMINIMIZER_H_
13 
14 #include <string>
15 
22 
23 namespace oops {
24 
25 /// GMRESR Minimizer
26 /*!
27  * Implements the GMRESR algorithm.
28  */
29 
30 // -----------------------------------------------------------------------------
31 
32 template<typename MODEL> class GMRESRMinimizer : public PrimalMinimizer<MODEL> {
37 
38  public:
39  const std::string classname() const override {return "GMRESRMinimizer";}
40  GMRESRMinimizer(const eckit::Configuration &, const CostFct_ & J): PrimalMinimizer<MODEL>(J) {}
42 
43  private:
44  double solve(CtrlInc_ &, const CtrlInc_ &,
45  const Hessian_ &, const Bmat_ &,
46  const int, const double) override;
47 };
48 
49 // =============================================================================
50 
51 template<typename MODEL>
53  const Hessian_ & hessian, const Bmat_ & B,
54  const int ninner, const double gnreduc) {
55 // Solve the linear system
56  double reduc = GMRESR(dx, rhs, hessian, B, ninner, gnreduc);
57  return reduc;
58 }
59 
60 // -----------------------------------------------------------------------------
61 
62 } // namespace oops
63 
64 #endif // OOPS_ASSIMILATION_GMRESRMINIMIZER_H_
double GMRESR(VECTOR &xx, const VECTOR &bb, const AMATRIX &A, const PMATRIX &precond, const int maxiter, const double tolerance)
Definition: GMRESR.h:63
BMatrix< MODEL > Bmat_
Cost Function.
Definition: CostFunction.h:56
The namespace for the main oops code.
Primal Minimizer.
CostFunction< MODEL > CostFct_
double solve(CtrlInc_ &, const CtrlInc_ &, const Hessian_ &, const Bmat_ &, const int, const double) override
GMRESR Minimizer.
The Hessian matrix: .
Definition: HessianMatrix.h:33
GMRESRMinimizer(const eckit::Configuration &, const CostFct_ &J)
HessianMatrix< MODEL > Hessian_
GMRESR solver for Ax=b.
const std::string classname() const override
The matrix.
Definition: BMatrix.h:27
ControlIncrement< MODEL > CtrlInc_