FV3 Bundle
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 <string>
15 #include <vector>
16 
17 #include <boost/scoped_ptr.hpp>
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> class SaddlePointMinimizer : public Minimizer<MODEL> {
48 // Eigen typedef SaddlePointLMPMatrix<MODEL> 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_(0) {}
55 
56  private:
57  CtrlInc_ * doMinimize(const eckit::Configuration &) override;
58 
59  const CostFct_ & J_;
60  boost::scoped_ptr<CtrlInc_> gradJb_;
61 // Eigen boost::scoped_ptr<LMP_> Pinv_;
62  std::vector< SaddlePointVector<MODEL> > xyVEC_;
63  std::vector< SaddlePointVector<MODEL> > pqVEC_;
64 };
65 
66 // =============================================================================
67 
68 template<typename MODEL>
70 SaddlePointMinimizer<MODEL>::doMinimize(const eckit::Configuration & config) {
71  int ninner = config.getInt("ninner");
72  int gnreduc = config.getDouble("gradient_norm_reduction");
73 
74 // if (gradJb_ == 0) 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> 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  boost::scoped_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> 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  Log::test() << "SaddlePointMinimizer: reduction in residual norm = " << reduc << std::endl;
121  Log::info() << "SaddlePointMinimizer output" << *dx << std::endl;
122 
123 // *gradJb_ = dx.lambda().dx();
124 // *gradJb_ *= -1.0;
125 
126  return dx;
127 }
128 
129 // -----------------------------------------------------------------------------
130 
131 } // namespace oops
132 
133 #endif // OOPS_ASSIMILATION_SADDLEPOINTMINIMIZER_H_
double GMRESR(VECTOR &xx, const VECTOR &bb, const AMATRIX &A, const PMATRIX &precond, const int maxiter, const double tolerance)
Definition: GMRESR.h:63
std::vector< SaddlePointVector< MODEL > > xyVEC_
FullGMRES solver for Ax=b.
Container of dual space vectors for all terms of the cost function.
Definition: DualVector.h:34
void dx(CtrlInc_ *dx)
Definition: DualVector.h:45
ControlIncrement< MODEL > CtrlInc_
CostFunction< MODEL > CostFct_
Cost Function.
Definition: CostFunction.h:56
program test
The namespace for the main oops code.
The preconditioner for the saddle-point minimizer.
The Saddle-point matrix.
SaddlePoint Minimizer.
subroutine, public info(self)
CtrlInc_ * doMinimize(const eckit::Configuration &) override
GMRESR solver for Ax=b.
void append(GeneralizedDepartures *)
Definition: DualVector.h:107
SaddlePointMinimizer(const eckit::Configuration &, const CostFct_ &J)
Control vector for the saddle point formulation.
A Minimizer knows how to minimize a cost function.
Definition: Minimizer.h:39
DualVector< MODEL > Multipliers_
const std::string classname() const override
std::vector< SaddlePointVector< MODEL > > pqVEC_
boost::scoped_ptr< CtrlInc_ > gradJb_