OOPS
SaddlePointMatrix.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_SADDLEPOINTMATRIX_H_
12 #define OOPS_ASSIMILATION_SADDLEPOINTMATRIX_H_
13 
14 #include <memory>
15 #include <utility>
16 
17 #include <boost/noncopyable.hpp>
18 
24 
25 namespace oops {
26 
27 /// The Saddle-point matrix.
28 /*!
29  * The solvers represent matrices as objects that implement a "multiply"
30  * method. This class defines objects that apply the saddle-point matrix.
31  */
32 
33 template<typename MODEL, typename OBS>
34 class SaddlePointMatrix : private boost::noncopyable {
38 
39  public:
40  explicit SaddlePointMatrix(const CostFct_ & j): j_(j) {}
41  void multiply(const SPVector_ &, SPVector_ &) const;
42 
43  private:
44  CostFct_ const & j_;
45 };
46 
47 // =============================================================================
48 
49 template<typename MODEL, typename OBS>
51  CtrlInc_ ww(j_.jb());
52 
53 // The three blocks below could be done in parallel
54 
55 // ADJ block
57  j_.zeroAD(ww);
58  z.dx(new CtrlInc_(j_.jb()));
59  j_.jb().initializeAD(z.dx(), x.lambda().dx(), costad);
60  for (unsigned jj = 0; jj < j_.nterms(); ++jj) {
61  j_.jterm(jj).computeCostAD(x.lambda().getv(jj), ww, costad);
62  }
63  j_.runADJ(ww, costad);
64  z.dx() += ww;
65  for (unsigned jj = 0; jj < j_.nterms(); ++jj) {
66  j_.jterm(jj).setPostProcAD();
67  }
68 
69 // TLM block
71  j_.jb().initializeTL(costtl);
72  for (unsigned jj = 0; jj < j_.nterms(); ++jj) {
73  j_.jterm(jj).setPostProcTL(x.dx(), costtl);
74  }
75  CtrlInc_ mdx(x.dx());
76  j_.runTLM(mdx, costtl);
77  z.lambda().clear();
78  z.lambda().dx(new CtrlInc_(j_.jb()));
79  j_.jb().finalizeTL(x.dx(), z.lambda().dx());
80  for (unsigned jj = 0; jj < j_.nterms(); ++jj) {
81  std::unique_ptr<GeneralizedDepartures> ztmp(j_.jterm(jj).newDualVector());
82  j_.jterm(jj).computeCostTL(x.dx(), *ztmp);
83  z.lambda().append(std::move(ztmp));
84  }
85 
86 // Diagonal block
88  diag.dx(new CtrlInc_(j_.jb()));
89  j_.jb().multiplyB(x.lambda().dx(), diag.dx());
90  for (unsigned jj = 0; jj < j_.nterms(); ++jj) {
91  diag.append(j_.jterm(jj).multiplyCovar(*x.lambda().getv(jj)));
92  }
93 
94 // The three blocks above could be done in parallel
95 
96  z.lambda() += diag;
97 }
98 
99 // -----------------------------------------------------------------------------
100 
101 } // namespace oops
102 
103 #endif // OOPS_ASSIMILATION_SADDLEPOINTMATRIX_H_
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
std::shared_ptr< const GeneralizedDepartures > getv(const unsigned) const
Definition: DualVector.h:129
Control model post processing.
void initializeTL(const Increment_ &dx, const util::DateTime &end, const util::Duration &step)
Tangent linear methods.
The Saddle-point matrix.
ControlIncrement< MODEL, OBS > CtrlInc_
void multiply(const SPVector_ &, SPVector_ &) const
CostFunction< MODEL, OBS > CostFct_
SaddlePointMatrix(const CostFct_ &j)
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.