OOPS
SaddlePointMinimizer.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_SADDLEPOINTMINIMIZER_H_
12 #define OOPS_ASSIMILATION_SADDLEPOINTMINIMIZER_H_
13 
14 #include <iomanip>
15 #include <memory>
16 #include <string>
17 #include <vector>
18 
19 #include "eckit/config/Configuration.h"
20 
27 // Eigen #include "oops/assimilation/SaddlePointLMPMatrix.h"
31 
32 #include "oops/util/Logger.h"
33 
34 namespace oops {
35 
36 /// SaddlePoint Minimizer
37 /*!
38  * Implements the SaddlePoint algorithm.
39  */
40 
41 // -----------------------------------------------------------------------------
42 
43 template<typename MODEL, typename OBS> class SaddlePointMinimizer : public Minimizer<MODEL, OBS> {
48 // Eigen typedef SaddlePointLMPMatrix<MODEL, OBS> LMP_;
49 
50  public:
51  const std::string classname() const override {return "SaddlePointMinimizer";}
52  SaddlePointMinimizer(const eckit::Configuration &, const CostFct_ & J)
53  : Minimizer_(J), J_(J), gradJb_() {}
55 
56  private:
57  CtrlInc_ * doMinimize(const eckit::Configuration &) override;
58 
59  const CostFct_ & J_;
60  std::unique_ptr<CtrlInc_> gradJb_;
61 // Eigen std::unique_ptr<LMP_> Pinv_;
62  std::vector< SaddlePointVector<MODEL, OBS> > xyVEC_;
63  std::vector< SaddlePointVector<MODEL, OBS> > pqVEC_;
64 };
65 
66 // =============================================================================
67 
68 template<typename MODEL, typename OBS>
70 SaddlePointMinimizer<MODEL, OBS>::doMinimize(const eckit::Configuration & config) {
71  int ninner = config.getInt("ninner");
72  double gnreduc = config.getDouble("gradient norm reduction");
73 
74 // if (!gradJb_) gradJb_.reset(new CtrlInc_(J_.jb()));
75 
76  Log::info() << "SaddlePointMinimizer: max iter = " << ninner
77  << ", requested norm reduction = " << gnreduc << std::endl;
78 
79 // Define saddle-point control vectors
80  Multipliers_ * pdxx = new Multipliers_();
81  pdxx->dx(new CtrlInc_(J_.jb()));
82  for (unsigned jj = 0; jj < J_.nterms(); ++jj) {
83  pdxx->append(J_.jterm(jj).newDualVector());
84  }
85  CtrlInc_ * tmp1 = new CtrlInc_(J_.jb());
86  SaddlePointVector<MODEL, OBS> spdx(tmp1, pdxx);
87 
88 // Compute RHS
89  Multipliers_ * pdfg = new Multipliers_();
90  CtrlInc_ * tmp3 = new CtrlInc_(J_.jb().getFirstGuess());
91  pdfg->dx(tmp3);
92  for (unsigned jj = 0; jj < J_.nterms(); ++jj) {
93  std::unique_ptr<GeneralizedDepartures> ww(J_.jterm(jj).newGradientFG());
94  pdfg->append(J_.jterm(jj).multiplyCovar(*ww));
95  }
96  CtrlInc_ * tmp2 = new CtrlInc_(J_.jb());
97  SaddlePointVector<MODEL, OBS> rhs(tmp2, pdfg);
98  rhs *= -1.0;
99 
100 // Define the matrices
102 
103 // Inexact constraint preconditioner
105 
106 // Initialize the limited memory preconditioner
107 // Eigen if (!Pinv_.get()) {
108 // Eigen Pinv_.reset(new LMP_(J_));
109 // Eigen }
110 
111 // Update the preconditioner
112 // Eigen Pinv_->setup(xyVEC_, pqVEC_);
113 
114 // Solve the linear system
115 // Eigen double reduc = FullGMRES(spdx, rhs, A, Pinv, ninner, gnreduc, pqVEC_, xyVEC_);
116  double reduc = GMRESR(spdx, rhs, A, Pinv, ninner, gnreduc);
117 
118  CtrlInc_ * dx = new CtrlInc_(spdx.dx());
119 
120  std::streamsize ss = Log::test().precision();
121  Log::test() << "SaddlePointMinimizer: reduction in residual norm = "
122  << std::setprecision(4) << reduc << std::setprecision(ss) << std::endl;
123  Log::info() << "SaddlePointMinimizer output" << *dx << std::endl;
124 
125 // *gradJb_ = dx.lambda().dx();
126 // *gradJb_ *= -1.0;
127 
128  return dx;
129 }
130 
131 // -----------------------------------------------------------------------------
132 
133 } // namespace oops
134 
135 #endif // OOPS_ASSIMILATION_SADDLEPOINTMINIMIZER_H_
GMRESR solver for Ax=b.
Cost Function.
Definition: CostFunction.h:53
Container of dual space vectors for all terms of the cost function.
Definition: DualVector.h:35
void dx(CtrlInc_ *dx)
Definition: DualVector.h:46
void append(std::unique_ptr< GeneralizedDepartures > &&)
Definition: DualVector.h:110
A Minimizer knows how to minimize a cost function.
Definition: Minimizer.h:37
The Saddle-point matrix.
SaddlePoint Minimizer.
DualVector< MODEL, OBS > Multipliers_
ControlIncrement< MODEL, OBS > CtrlInc_
CostFunction< MODEL, OBS > CostFct_
std::vector< SaddlePointVector< MODEL, OBS > > pqVEC_
SaddlePointMinimizer(const eckit::Configuration &, const CostFct_ &J)
std::vector< SaddlePointVector< MODEL, OBS > > xyVEC_
Minimizer< MODEL, OBS > Minimizer_
CtrlInc_ * doMinimize(const eckit::Configuration &) override
std::unique_ptr< CtrlInc_ > gradJb_
const std::string classname() const override
The preconditioner for the saddle-point minimizer.
Control vector for the saddle point formulation.
const CtrlInc_ & dx() const
Accessor method to get the dx_ component.
The namespace for the main oops code.
double GMRESR(VECTOR &xx, const VECTOR &bb, const AMATRIX &A, const PMATRIX &precond, const int maxiter, const double tolerance)
Definition: GMRESR.h:64
FullGMRES solver for Ax=b.