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