FV3 Bundle
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 <boost/noncopyable.hpp>
15 
22 
23 namespace oops {
24  template<typename MODEL> class JqTermTLAD;
25 
26 /// The Saddle-point matrix.
27 /*!
28  * The solvers represent matrices as objects that implement a "multiply"
29  * method. This class defines objects that apply the saddle-point matrix.
30  */
31 
32 template<typename MODEL>
33 class SaddlePointMatrix : private boost::noncopyable {
39 
40  public:
41  explicit SaddlePointMatrix(const CostFct_ & j): j_(j) {}
42  void multiply(const SPVector_ &, SPVector_ &) const;
43 
44  private:
45  CostFct_ const & j_;
46 };
47 
48 // =============================================================================
49 
50 template<typename MODEL>
52  SPVector_ & z) const {
53  CtrlInc_ ww(j_.jb());
54 
55 // The three blocks below could be done in parallel
56 
57 // ADJ block
59  j_.zeroAD(ww);
60  z.dx(new CtrlInc_(j_.jb()));
61  JqTermTLAD_ * jqad = j_.jb().initializeAD(z.dx(), x.lambda().dx());
62  costad.enrollProcessor(jqad);
63  for (unsigned jj = 0; jj < j_.nterms(); ++jj) {
64  costad.enrollProcessor(j_.jterm(jj).setupAD(x.lambda().getv(jj), ww));
65  }
66  j_.runADJ(ww, costad);
67  z.dx() += ww;
68 
69 // TLM block
71  JqTermTLAD_ * jqtl = j_.jb().initializeTL();
72  costtl.enrollProcessor(jqtl);
73  for (unsigned jj = 0; jj < j_.nterms(); ++jj) {
74  costtl.enrollProcessor(j_.jterm(jj).setupTL(x.dx()));
75  }
76  CtrlInc_ mdx(x.dx());
77  j_.runTLM(mdx, costtl);
78  z.lambda().clear();
79  z.lambda().dx(new CtrlInc_(j_.jb()));
80  j_.jb().finalizeTL(jqtl, x.dx(), z.lambda().dx());
81  for (unsigned jj = 0; jj < j_.nterms(); ++jj) {
82  z.lambda().append(costtl.releaseOutputFromTL(jj+1));
83  }
84 
85 // Diagonal block
86  DualVector<MODEL> diag;
87  diag.dx(new CtrlInc_(j_.jb()));
88  j_.jb().multiplyB(x.lambda().dx(), diag.dx());
89  for (unsigned jj = 0; jj < j_.nterms(); ++jj) {
90  diag.append(j_.jterm(jj).multiplyCovar(*x.lambda().getv(jj)));
91  }
92 
93 // The three blocks above could be done in parallel
94 
95  z.lambda() += diag;
96 }
97 
98 // -----------------------------------------------------------------------------
99 
100 } // namespace oops
101 
102 #endif // OOPS_ASSIMILATION_SADDLEPOINTMATRIX_H_
GeneralizedDepartures * releaseOutputFromTL(unsigned int ii)
Get TL dual space output.
Container of dual space vectors for all terms of the cost function.
Definition: DualVector.h:34
SaddlePointMatrix(const CostFct_ &j)
void dx(CtrlInc_ *dx)
Definition: DualVector.h:45
void initializeTL(const Increment_ &dx, const util::DateTime &end, const util::Duration &step)
Tangent linear methods.
Definition: PostBaseTLAD.h:66
Cost Function.
Definition: CostFunction.h:56
l_size ! loop over number of fields ke do j
The namespace for the main oops code.
SaddlePointVector< MODEL > SPVector_
ControlIncrement< MODEL > CtrlInc_
The Saddle-point matrix.
void initializeAD(Increment_ &dx, const util::DateTime &bgn, const util::Duration &step)
Adjoint methods.
Definition: PostBaseTLAD.h:84
Control model post processing.
Increment< MODEL > Increment_
void enrollProcessor(PostBaseTLAD_ *pp)
void multiply(const SPVector_ &, SPVector_ &) const
JqTermTLAD< MODEL > JqTermTLAD_
Increment Class: Difference between two states.
const Multipliers_ & lambda() const
Accessor method to get the lambda_ component.
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.
CostFunction< MODEL > CostFct_
boost::shared_ptr< const GeneralizedDepartures > getv(const unsigned) const
Definition: DualVector.h:127