OOPS
SaddlePointPrecondMatrix.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_SADDLEPOINTPRECONDMATRIX_H_
12 #define OOPS_ASSIMILATION_SADDLEPOINTPRECONDMATRIX_H_
13 
14 #include <boost/noncopyable.hpp>
15 
21 
22 namespace oops {
23 
24  /// The preconditioner for the saddle-point minimizer.
25  /*!
26  * The solvers represent matrices as objects that implement a "multiply"
27  * method. This class defines objects that apply the saddle-point matrix.
28  */
29 
30 // -----------------------------------------------------------------------------
31 template<typename MODEL, typename OBS>
32 class SaddlePointPrecondMatrix : private boost::noncopyable {
37 
38  public:
39  explicit SaddlePointPrecondMatrix(const CostFct_ & j);
40  void multiply(const SPVector_ &, SPVector_ &) const;
41 
42  private:
43  void Lhatinv(const CtrlInc_ &, CtrlInc_ &, const int) const;
44  void Lhatinvt(const CtrlInc_ &, CtrlInc_ &, const int) const;
45 
46  const CostFctWeak_ & j_;
47  const bool idmodel_;
48 };
49 
50 // =============================================================================
51 
52 template<typename MODEL, typename OBS>
54  : j_(dynamic_cast<const CostFctWeak_ &>(j)),
55  idmodel_(false)
56 {}
57 
58 // -----------------------------------------------------------------------------
59 
60 template<typename MODEL, typename OBS>
62  z.lambda().clear();
63  z.lambda().dx(new CtrlInc_(j_.jb()));
64  z.dx(new CtrlInc_(j_.jb()));
65 
66  int norder = 1; // truncation order for Lhatinv and Lhatinvt
67 
68  Lhatinvt(x.dx(), z.lambda().dx(), norder);
69 
70  CtrlInc_ l(z.lambda().dx());
71  j_.jb().multiplyB(z.lambda().dx(), l);
72  l *= -1.0;
73  l += x.lambda().dx();
74 
75  Lhatinv(l, z.dx(), norder);
76 
77  for (unsigned jj = 0; jj < j_.nterms(); ++jj) {
78  z.lambda().append(j_.jterm(jj).multiplyCoInv(*x.lambda().getv(jj)));
79  }
80 }
81 
82 // -----------------------------------------------------------------------------
83 
84 template<typename MODEL, typename OBS>
86  const int norder) const {
87 // Approximate L=I + (I-L) + (I-L)^2 + ... + (I-L)^norder
88  CtrlInc_ ww(xx);
89 
90  zz = xx;
91  for (int jj = 0; jj < norder; ++jj) {
92  j_.runTLM(ww, idmodel_);
93  ww.shift_forward();
94  zz += ww;
95  }
96 }
97 
98 // -----------------------------------------------------------------------------
99 
100 template<typename MODEL, typename OBS>
102  const int norder) const {
103 // Approximate L'=I + (I-L') + (I-L')^2 + ... + (I-L')^norder
104  CtrlInc_ ww(xx);
105 
106  zz = xx;
107  for (int jj = 0; jj < norder; ++jj) {
108  ww.shift_backward();
109  j_.runADJ(ww, idmodel_);
110  zz += ww;
111  }
112 }
113 
114 // -----------------------------------------------------------------------------
115 } // namespace oops
116 
117 #endif // OOPS_ASSIMILATION_SADDLEPOINTPRECONDMATRIX_H_
Weak Constraint 4D-Var Cost Function.
Definition: CostFctWeak.h:53
Cost Function.
Definition: CostFunction.h:53
void dx(CtrlInc_ *dx)
Definition: DualVector.h:46
void append(std::unique_ptr< GeneralizedDepartures > &&)
Definition: DualVector.h:110
std::shared_ptr< const GeneralizedDepartures > getv(const unsigned) const
Definition: DualVector.h:129
The preconditioner for the saddle-point minimizer.
void Lhatinv(const CtrlInc_ &, CtrlInc_ &, const int) const
void Lhatinvt(const CtrlInc_ &, CtrlInc_ &, const int) const
ControlIncrement< MODEL, OBS > CtrlInc_
CostFctWeak< MODEL, OBS > CostFctWeak_
void multiply(const SPVector_ &, SPVector_ &) const
CostFunction< MODEL, OBS > CostFct_
SaddlePointVector< MODEL, OBS > SPVector_
Control vector for the saddle point formulation.
const Multipliers_ & lambda() const
Accessor method to get the lambda_ component.
const CtrlInc_ & dx() const
Accessor method to get the dx_ component.
The namespace for the main oops code.