IODA Bundle
ObsErrorDiag.h
Go to the documentation of this file.
1 /*
2  * (C) Copyright 2009-2016 ECMWF.
3  * (C) Copyright 2021 UCAR.
4  *
5  * This software is licensed under the terms of the Apache Licence Version 2.0
6  * which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
7  * In applying this licence, ECMWF does not waive the privileges and immunities
8  * granted to it by virtue of its status as an intergovernmental organisation nor
9  * does it submit to any jurisdiction.
10  */
11 
12 #ifndef OOPS_GENERIC_OBSERRORDIAG_H_
13 #define OOPS_GENERIC_OBSERRORDIAG_H_
14 
15 #include <sstream>
16 #include <string>
17 
18 #include "eckit/config/Configuration.h"
22 #include "oops/util/Logger.h"
23 #include "oops/util/parameters/Parameter.h"
24 #include "oops/util/parameters/Parameters.h"
25 
26 namespace oops {
27 
28 /// \brief Parameters for diagonal obs errors
29 class ObsErrorDiagParameters : public Parameters {
30  OOPS_CONCRETE_PARAMETERS(ObsErrorDiagParameters, Parameters)
31  public:
32  /// perturbation amplitude multiplier
33  Parameter<double> pert{"random amplitude", 1.0, this};
34 };
35 
36 // -----------------------------------------------------------------------------
37 /// \brief Diagonal observation error covariance matrix.
38 template<typename OBS>
39 class ObsErrorDiag : public ObsErrorBase<OBS> {
42 
43  public:
44  ObsErrorDiag(const eckit::Configuration &, const ObsSpace_ &);
45 
46 /// Update after obs errors potentially changed
47  void update(const ObsVector_ &) override;
48 
49 /// Multiply a Departure by \f$R\f$
50  void multiply(ObsVector_ &) const override;
51 
52 /// Multiply a Departure by \f$R^{-1}\f$
53  void inverseMultiply(ObsVector_ &) const override;
54 
55 /// Generate random perturbation
56  void randomize(ObsVector_ &) const override;
57 
58 /// Save obs errors
59  void save(const std::string &) const override;
60 
61 /// Get mean std deviation of errors for Jo table
62  double getRMSE() const override {return stddev_.rms();}
63 
64 /// Get obs errors std deviation
65  ObsVector_ obserrors() const override {return stddev_;}
66 
67 /// Get inverseVariance
68  ObsVector_ inverseVariance() const override {return inverseVariance_;}
69 
70  private:
71  void print(std::ostream &) const override;
75 };
76 
77 // =============================================================================
78 
79 template<typename OBS>
80 ObsErrorDiag<OBS>::ObsErrorDiag(const eckit::Configuration & conf, const ObsSpace_ & obsgeom)
81  : stddev_(obsgeom, "ObsError"), inverseVariance_(obsgeom)
82 {
83  options_.deserialize(conf);
87  Log::trace() << "ObsErrorDiag:ObsErrorDiag constructed nobs = " << stddev_.nobs() << std::endl;
88 }
89 
90 // -----------------------------------------------------------------------------
91 
92 template<typename OBS>
93 void ObsErrorDiag<OBS>::update(const ObsVector_ & obserr) {
94  stddev_ = obserr;
95  inverseVariance_ = stddev_;
96  inverseVariance_ *= stddev_;
97  inverseVariance_.invert();
98  Log::info() << "ObsErrorDiag covariance updated " << stddev_.nobs() << std::endl;
99 }
100 
101 // -----------------------------------------------------------------------------
102 
103 template<typename OBS>
105  dy /= inverseVariance_;
106 }
107 
108 // -----------------------------------------------------------------------------
109 
110 template<typename OBS>
112  dy *= inverseVariance_;
113 }
114 
115 // -----------------------------------------------------------------------------
116 
117 template<typename OBS>
119  dy.random();
120  dy *= stddev_;
121  dy *= options_.pert;
122 }
123 
124 // -----------------------------------------------------------------------------
125 
126 template<typename OBS>
127 void ObsErrorDiag<OBS>::save(const std::string & name) const {
128  stddev_.save(name);
129 }
130 
131 // -----------------------------------------------------------------------------
132 
133 template<typename OBS>
134 void ObsErrorDiag<OBS>::print(std::ostream & os) const {
135  os << "Diagonal observation error covariance" << std::endl << stddev_;
136 }
137 
138 // -----------------------------------------------------------------------------
139 
140 
141 } // namespace oops
142 
143 #endif // OOPS_GENERIC_OBSERRORDIAG_H_
Base class for generic implementations of observation error covariance matrices.
Diagonal observation error covariance matrix.
Definition: ObsErrorDiag.h:39
double getRMSE() const override
Get mean std deviation of errors for Jo table.
Definition: ObsErrorDiag.h:62
void save(const std::string &) const override
Save obs errors.
Definition: ObsErrorDiag.h:127
void print(std::ostream &) const override
Definition: ObsErrorDiag.h:134
ObsVector_ inverseVariance() const override
Get inverseVariance.
Definition: ObsErrorDiag.h:68
void inverseMultiply(ObsVector_ &) const override
Multiply a Departure by .
Definition: ObsErrorDiag.h:111
void multiply(ObsVector_ &) const override
Multiply a Departure by .
Definition: ObsErrorDiag.h:104
ObsVector_ stddev_
Definition: ObsErrorDiag.h:72
ObsVector_ inverseVariance_
Definition: ObsErrorDiag.h:73
ObsSpace< OBS > ObsSpace_
Definition: ObsErrorDiag.h:40
void randomize(ObsVector_ &) const override
Generate random perturbation.
Definition: ObsErrorDiag.h:118
ObsVector_ obserrors() const override
Get obs errors std deviation.
Definition: ObsErrorDiag.h:65
ObsErrorDiag(const eckit::Configuration &, const ObsSpace_ &)
Definition: ObsErrorDiag.h:80
ObsErrorDiagParameters options_
Definition: ObsErrorDiag.h:74
void update(const ObsVector_ &) override
Update after obs errors potentially changed.
Definition: ObsErrorDiag.h:93
ObsVector< OBS > ObsVector_
Definition: ObsErrorDiag.h:41
Parameters for diagonal obs errors.
Definition: ObsErrorDiag.h:29
Parameter< double > pert
perturbation amplitude multiplier
Definition: ObsErrorDiag.h:33
unsigned int nobs() const
number of non-masked out observations (across all MPI tasks)
The namespace for the main oops code.