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 <boost/noncopyable.hpp>
15
16
#include "
oops/assimilation/ControlIncrement.h
"
17
#include "
oops/assimilation/CostFunction.h
"
18
#include "
oops/assimilation/DualVector.h
"
19
#include "
oops/assimilation/SaddlePointVector.h
"
20
#include "
oops/base/PostProcessorTLAD.h
"
21
22
namespace
oops
{
23
template
<
typename
MODEL>
class
JqTermTLAD;
24
25
/// The Saddle-point matrix.
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
template
<
typename
MODEL,
typename
OBS>
32
class
SaddlePointMatrix
:
private
boost::noncopyable {
33
typedef
ControlIncrement<MODEL, OBS>
CtrlInc_
;
34
typedef
CostFunction<MODEL, OBS>
CostFct_
;
35
typedef
SaddlePointVector<MODEL, OBS>
SPVector_
;
36
typedef
JqTermTLAD<MODEL>
JqTermTLAD_
;
37
38
public
:
39
explicit
SaddlePointMatrix
(
const
CostFct_
& j):
j_
(j) {}
40
void
multiply
(
const
SPVector_
&,
SPVector_
&)
const
;
41
42
private
:
43
CostFct_
const
&
j_
;
44
};
45
46
// =============================================================================
47
48
template
<
typename
MODEL,
typename
OBS>
49
void
SaddlePointMatrix<MODEL, OBS>::multiply
(
const
SPVector_
& x,
50
SPVector_
& z)
const
{
51
CtrlInc_
ww(j_.jb());
52
53
// The three blocks below could be done in parallel
54
55
// ADJ block
56
PostProcessorTLAD<MODEL>
costad;
57
j_.zeroAD(ww);
58
z.
dx
(
new
CtrlInc_
(j_.jb()));
59
JqTermTLAD_
* jqad = j_.jb().initializeAD(z.
dx
(), x.
lambda
().
dx
());
60
costad.
enrollProcessor
(jqad);
61
for
(
unsigned
jj = 0; jj < j_.nterms(); ++jj) {
62
costad.
enrollProcessor
(j_.jterm(jj).setupAD(x.
lambda
().
getv
(jj), ww));
63
}
64
j_.runADJ(ww, costad);
65
z.
dx
() += ww;
66
67
// TLM block
68
PostProcessorTLAD<MODEL>
costtl;
69
JqTermTLAD_
* jqtl = j_.jb().initializeTL();
70
costtl.
enrollProcessor
(jqtl);
71
for
(
unsigned
jj = 0; jj < j_.nterms(); ++jj) {
72
costtl.
enrollProcessor
(j_.jterm(jj).setupTL(x.
dx
()));
73
}
74
CtrlInc_
mdx(x.
dx
());
75
j_.runTLM(mdx, costtl);
76
z.
lambda
().
clear
();
77
z.
lambda
().
dx
(
new
CtrlInc_
(j_.jb()));
78
j_.jb().finalizeTL(jqtl, x.
dx
(), z.
lambda
().
dx
());
79
for
(
unsigned
jj = 0; jj < j_.nterms(); ++jj) {
80
z.
lambda
().
append
(costtl.
releaseOutputFromTL
(jj+1));
81
}
82
83
// Diagonal block
84
DualVector<MODEL, OBS>
diag;
85
diag.
dx
(
new
CtrlInc_
(j_.jb()));
86
j_.jb().multiplyB(x.
lambda
().
dx
(), diag.
dx
());
87
for
(
unsigned
jj = 0; jj < j_.nterms(); ++jj) {
88
diag.
append
(j_.jterm(jj).multiplyCovar(*x.
lambda
().
getv
(jj)));
89
}
90
91
// The three blocks above could be done in parallel
92
93
z.
lambda
() += diag;
94
}
95
96
// -----------------------------------------------------------------------------
97
98
}
// namespace oops
99
100
#endif // OOPS_ASSIMILATION_SADDLEPOINTMATRIX_H_
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::SaddlePointMatrix::j_
CostFct_ const & j_
Definition:
SaddlePointMatrix.h:43
oops::SaddlePointMatrix::CtrlInc_
ControlIncrement< MODEL, OBS > CtrlInc_
Definition:
SaddlePointMatrix.h:33
SaddlePointVector.h
CostFunction.h
oops::SaddlePointMatrix::multiply
void multiply(const SPVector_ &, SPVector_ &) const
Definition:
SaddlePointMatrix.h:49
oops::DualVector
Container of dual space vectors for all terms of the cost function.
Definition:
DualVector.h:34
oops::DualVector::dx
void dx(CtrlInc_ *dx)
Definition:
DualVector.h:45
oops::SaddlePointMatrix
The Saddle-point matrix.
Definition:
SaddlePointMatrix.h:32
oops::ControlIncrement
Definition:
ControlIncrement.h:46
oops::DualVector::append
void append(std::unique_ptr< GeneralizedDepartures > &&)
Definition:
DualVector.h:107
oops::JqTermTLAD
Definition:
CostJb3D.h:29
oops::SaddlePointVector::dx
const CtrlInc_ & dx() const
Accessor method to get the dx_ component.
Definition:
SaddlePointVector.h:46
oops::PostProcessorTLAD::releaseOutputFromTL
std::unique_ptr< GeneralizedDepartures > releaseOutputFromTL(unsigned int ii)
Get TL dual space output.
Definition:
PostProcessorTLAD.h:95
oops::SaddlePointMatrix::SaddlePointMatrix
SaddlePointMatrix(const CostFct_ &j)
Definition:
SaddlePointMatrix.h:39
oops::PostProcessorTLAD::enrollProcessor
void enrollProcessor(PostBaseTLAD_ *pp)
Definition:
PostProcessorTLAD.h:43
oops::SaddlePointMatrix::CostFct_
CostFunction< MODEL, OBS > CostFct_
Definition:
SaddlePointMatrix.h:34
oops::SaddlePointMatrix::SPVector_
SaddlePointVector< MODEL, OBS > SPVector_
Definition:
SaddlePointMatrix.h:35
oops::SaddlePointVector
Control vector for the saddle point formulation.
Definition:
SaddlePointVector.h:30
oops::PostProcessorTLAD
Control model post processing.
Definition:
PostProcessorTLAD.h:33
oops::SaddlePointMatrix::JqTermTLAD_
JqTermTLAD< MODEL > JqTermTLAD_
Definition:
SaddlePointMatrix.h:36
oops::DualVector::getv
std::shared_ptr< const GeneralizedDepartures > getv(const unsigned) const
Definition:
DualVector.h:126
ControlIncrement.h
oops::CostFunction
Cost Function.
Definition:
CostFunction.h:53
PostProcessorTLAD.h
DualVector.h
fv3-bundle
oops
src
oops
assimilation
SaddlePointMatrix.h
Generated on Sun Oct 25 2020 12:42:59 for OOPS by
1.8.18