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