OOPS
PCGMinimizer.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_PCGMINIMIZER_H_
12 #define OOPS_ASSIMILATION_PCGMINIMIZER_H_
13 
14 #include <string>
15 
20 #include "oops/assimilation/PCG.h"
22 
23 namespace oops {
24 
25 /// PCG Minimizer
26 /*!
27  * Implements the standard Preconditioned Conjugate Gradients algorithm.
28  */
29 
30 // -----------------------------------------------------------------------------
31 
32 template<typename MODEL, typename OBS> class PCGMinimizer : public PrimalMinimizer<MODEL, OBS> {
37 
38  public:
39  const std::string classname() const override {return "PCGMinimizer";}
40  PCGMinimizer(const eckit::Configuration &, const CostFct_ & J): PrimalMinimizer<MODEL, OBS>(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, typename OBS>
53  const Hessian_ & hessian, const Bmat_ & B,
54  const int ninner, const double gnreduc) {
55 // Solve the linear system
56  double reduc = PCG(dx, rhs, hessian, B, ninner, gnreduc);
57  return reduc;
58 }
59 
60 // -----------------------------------------------------------------------------
61 
62 } // namespace oops
63 
64 #endif // OOPS_ASSIMILATION_PCGMINIMIZER_H_
oops::PCGMinimizer::~PCGMinimizer
~PCGMinimizer()
Definition: PCGMinimizer.h:41
oops::PrimalMinimizer
Primal Minimizer.
Definition: PrimalMinimizer.h:33
oops
The namespace for the main oops code.
Definition: ErrorCovarianceL95.cc:22
PCG.h
Preconditioned Conjugate Gradients solver.
HessianMatrix.h
oops::PCGMinimizer::solve
double solve(CtrlInc_ &, const CtrlInc_ &, const Hessian_ &, const Bmat_ &, const int, const double) override
Definition: PCGMinimizer.h:52
CostFunction.h
PrimalMinimizer.h
oops::PCGMinimizer::PCGMinimizer
PCGMinimizer(const eckit::Configuration &, const CostFct_ &J)
Definition: PCGMinimizer.h:40
oops::ControlIncrement
Definition: ControlIncrement.h:46
oops::PCGMinimizer::Bmat_
BMatrix< MODEL, OBS > Bmat_
Definition: PCGMinimizer.h:33
oops::PCGMinimizer::classname
const std::string classname() const override
Definition: PCGMinimizer.h:39
oops::PCGMinimizer::CtrlInc_
ControlIncrement< MODEL, OBS > CtrlInc_
Definition: PCGMinimizer.h:35
oops::PCGMinimizer
PCG Minimizer.
Definition: PCGMinimizer.h:32
oops::HessianMatrix
The Hessian matrix: .
Definition: HessianMatrix.h:31
BMatrix.h
oops::PCG
double PCG(VECTOR &x, const VECTOR &b, const AMATRIX &A, const PMATRIX &precond, const int maxiter, const double tolerance)
Definition: PCG.h:63
ControlIncrement.h
oops::PCGMinimizer::Hessian_
HessianMatrix< MODEL, OBS > Hessian_
Definition: PCGMinimizer.h:36
oops::BMatrix
The matrix.
Definition: BMatrix.h:27
oops::CostFunction
Cost Function.
Definition: CostFunction.h:53
oops::PCGMinimizer::CostFct_
CostFunction< MODEL, OBS > CostFct_
Definition: PCGMinimizer.h:34