UFO
ObsErrorCrossVarCov.h
Go to the documentation of this file.
1 /*
2  * (C) Copyright 2021 UCAR.
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  */
7 
8 #ifndef UFO_ERRORS_OBSERRORCROSSVARCOV_H_
9 #define UFO_ERRORS_OBSERRORCROSSVARCOV_H_
10 
11 #include <Eigen/Core>
12 
13 #include <memory>
14 #include <string>
15 
16 #include "ioda/ObsVector.h"
17 
18 #include "oops/interface/ObsErrorBase.h"
19 #include "oops/util/parameters/Parameters.h"
20 #include "oops/util/parameters/RequiredParameter.h"
21 
22 #include "ufo/ObsTraits.h"
23 
24 namespace eckit {
25  class Configuration;
26 }
27 
28 namespace ioda {
29  class ObsSpace;
30 }
31 
32 namespace ufo {
33 
34 /// \brief Parameters for obs errors with cross-variable correlations
35 class ObsErrorCrossVarCovParameters : public oops::ObsErrorParametersBase {
36  OOPS_CONCRETE_PARAMETERS(ObsErrorCrossVarCovParameters, ObsErrorParametersBase)
37  public:
38  /// Input file containing correlations
39  oops::RequiredParameter<std::string> inputFile{"input file", this};
40 };
41 
42 // -----------------------------------------------------------------------------
43 /// \brief Observation error covariance matrix with cross-variable
44 /// correlations.
45 /// \details Correlations are the same at all locations and are read
46 /// from the file specified in the configuration. Obs error standard
47 /// deviations are read from ObsSpace as ObsError group.
48 /// Full observation error covariance matrix is R = D^{1/2} * C * D^{1/2}
49 /// where D^{1/2} is a diagonal matrix with stddev_ (ObsError group)
50 /// on the diagonal, and C is the correlation matrix.
51 class ObsErrorCrossVarCov : public oops::interface::ObsErrorBase<ObsTraits> {
52  public:
53  /// The type of parameters passed to the constructor.
54  /// This typedef is used by the ObsErrorFactory.
56 
57  /// Initialize observation errors
58  ObsErrorCrossVarCov(const Parameters_ &, ioda::ObsSpace &,
59  const eckit::mpi::Comm &timeComm);
60 
61  /// Update obs error standard deviations to be equal to \p stddev
62  void update(const ioda::ObsVector & stddev) override;
63 
64  /// Multiply \p y by this observation error covariance
65  /// Computed as R * dy = D^{1/2} * C * D^{1/2} * dy
66  /// where D^{1/2} - diagonal matrix with stddev_ on the diagonal
67  /// C - correlations
68  void multiply(ioda::ObsVector & y) const override;
69 
70  /// Multiply \p y by inverse of this observation error covariance
71  /// Computed as R^{-1} * dy = D^{-1/2} * C^{-1] * D^{-1/2} * dy
72  /// where D^{1/2} - diagonal matrix with stddev_ on the diagonal
73  /// C - correlations
74  void inverseMultiply(ioda::ObsVector & y) const override;
75 
76  /// Generate \p y as a random perturbation
77  void randomize(ioda::ObsVector & y) const override;
78 
79  /// Save obs error standard deviations under \p name group name
80  void save(const std::string & name) const override;
81 
82  /// Return RMS of obs error standard deviations
83  double getRMSE() const override {return stddev_.rms();}
84 
85  /// Return obs errors std deviation
86  std::unique_ptr<ioda::ObsVector> getObsErrors() const override;
87 
88  /// Return inverse of obs error variance
89  std::unique_ptr<ioda::ObsVector> getInverseVariance() const override;
90 
91  private:
92  /// Print covariance details (for logging)
93  void print(std::ostream &) const override;
94  /// Observation error standard deviations
95  ioda::ObsVector stddev_;
96  /// Correlations between variables
97  Eigen::MatrixXd varcorrelations_;
98 };
99 
100 // -----------------------------------------------------------------------------
101 
102 } // namespace ufo
103 
104 #endif // UFO_ERRORS_OBSERRORCROSSVARCOV_H_
Observation error covariance matrix with cross-variable correlations.
ObsErrorCrossVarCov(const Parameters_ &, ioda::ObsSpace &, const eckit::mpi::Comm &timeComm)
Initialize observation errors.
void randomize(ioda::ObsVector &y) const override
Generate y as a random perturbation.
void multiply(ioda::ObsVector &y) const override
std::unique_ptr< ioda::ObsVector > getObsErrors() const override
Return obs errors std deviation.
ObsErrorCrossVarCovParameters Parameters_
Eigen::MatrixXd varcorrelations_
Correlations between variables.
void save(const std::string &name) const override
Save obs error standard deviations under name group name.
void update(const ioda::ObsVector &stddev) override
Update obs error standard deviations to be equal to stddev.
ioda::ObsVector stddev_
Observation error standard deviations.
void inverseMultiply(ioda::ObsVector &y) const override
void print(std::ostream &) const override
Print covariance details (for logging)
std::unique_ptr< ioda::ObsVector > getInverseVariance() const override
Return inverse of obs error variance.
double getRMSE() const override
Return RMS of obs error standard deviations.
Parameters for obs errors with cross-variable correlations.
oops::RequiredParameter< std::string > inputFile
Input file containing correlations.
Forward declarations.
Definition: ObsAodExt.h:21
Forward declarations.
Definition: ObsAodExt.h:25
Definition: RunCRTM.h:27