FV3 Bundle
IPCGMinimizer.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_IPCGMINIMIZER_H_
12 #define OOPS_ASSIMILATION_IPCGMINIMIZER_H_
13 
14 #include <string>
15 
20 #include "oops/assimilation/IPCG.h"
22 
23 namespace oops {
24 
25 /// IPCG Minimizer
26 /*!
27  * Implements the Golub-Ye Inexact-Preconditioned Conjugate Gradients algorithm.
28  */
29 
30 // -----------------------------------------------------------------------------
31 
32 template<typename MODEL> class IPCGMinimizer : public PrimalMinimizer<MODEL> {
37 
38  public:
39  const std::string classname() const override {return "IPCGMinimizer";}
40  IPCGMinimizer(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 = IPCG(dx, rhs, hessian, B, ninner, gnreduc);
57  reduc = round(100000.0 * reduc)/100000.0; // Reducing precision for test, NOT GOOD!!!
58  return reduc;
59 }
60 
61 // -----------------------------------------------------------------------------
62 
63 } // namespace oops
64 
65 #endif // OOPS_ASSIMILATION_IPCGMINIMIZER_H_
HessianMatrix< MODEL > Hessian_
Definition: IPCGMinimizer.h:36
double IPCG(VECTOR &x, const VECTOR &b, const AMATRIX &A, const PMATRIX &precond, const int maxiter, const double tolerance)
Definition: IPCG.h:62
Cost Function.
Definition: CostFunction.h:56
IPCGMinimizer(const eckit::Configuration &, const CostFct_ &J)
Definition: IPCGMinimizer.h:40
The namespace for the main oops code.
Inexact-Preconditioned Conjugate Gradients solver.
Primal Minimizer.
BMatrix< MODEL > Bmat_
Definition: IPCGMinimizer.h:33
double solve(CtrlInc_ &, const CtrlInc_ &, const Hessian_ &, const Bmat_ &, const int, const double) override
Definition: IPCGMinimizer.h:52
ControlIncrement< MODEL > CtrlInc_
Definition: IPCGMinimizer.h:35
The Hessian matrix: .
Definition: HessianMatrix.h:33
The matrix.
Definition: BMatrix.h:27
const std::string classname() const override
Definition: IPCGMinimizer.h:39
IPCG Minimizer.
Definition: IPCGMinimizer.h:32
CostFunction< MODEL > CostFct_
Definition: IPCGMinimizer.h:34