UFO
src/ufo/ObsBias.h
Go to the documentation of this file.
1 /*
2  * (C) Copyright 2017-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_OBSBIAS_H_
9 #define UFO_OBSBIAS_H_
10 
11 #include <Eigen/Core>
12 
13 #include <memory>
14 #include <string>
15 #include <vector>
16 
17 #include "eckit/config/LocalConfiguration.h"
18 
19 #include "oops/base/Variables.h"
20 #include "oops/util/ObjectCounter.h"
21 #include "oops/util/Printable.h"
22 
23 #include "ufo/ObsBiasParameters.h"
25 
26 namespace oops {
27  class Variables;
28 }
29 
30 namespace ioda {
31  class ObsSpace;
32  class ObsVector;
33 }
34 
35 namespace ufo {
36  class GeoVals;
37  class ObsBiasIncrement;
38  class ObsDiagnostics;
39 
40 /// Class to handle observation bias correction coefficients
41 /// \details contains information on what predictors are used for bias
42 /// correction application
43 class ObsBias : public util::Printable,
44  private util::ObjectCounter<ObsBias> {
45  public:
47 
48  static const std::string classname() {return "ufo::ObsBias";}
49 
50  ObsBias(ioda::ObsSpace &, const Parameters_ &);
51  ObsBias(const ObsBias &, const bool);
52 
54  ObsBias & operator=(const ObsBias &);
55 
56  /// Read bias correction coefficients from file
57  void read(const Parameters_ &);
58  void write(const Parameters_ &) const;
59  double norm() const;
60  std::size_t size() const {return biascoeffs_.size();}
61 
62  /// Return the coefficient of predictor \p jpred for variable \p jvar.
63  ///
64  /// Note: \p jpred may be the index of a static or a variable predictor.
65  double operator()(size_t jpred, size_t jvar) const {
66  return jpred < numStaticPredictors_ ?
67  1.0 : biascoeffs_[index(jpred - numStaticPredictors_, jvar)];
68  }
69 
70  /// Return bias correction coefficients (for *variable* predictors)
71  const Eigen::VectorXd & data() const {return biascoeffs_;}
72 
73  // Required variables
74  const oops::Variables & requiredVars() const {return geovars_;}
75  const oops::Variables & requiredHdiagnostics() const {return hdiags_;}
76  const std::vector<std::string> & requiredPredictors() const {return prednames_;}
77 
78  /// Return a reference to the vector of all (static and variable) predictors.
79  const Predictors & predictors() const {return predictors_;}
80 
81  /// Return the vector of variable predictors.
82  std::vector<std::shared_ptr<const PredictorBase>> variablePredictors() const;
83 
84  /// Return the list of bias-corrected variables.
85  const oops::Variables & correctedVars() const {return vars_;}
86 
87  /// Set all variable predictors coefficients to zero (used in the test)
88  void zero();
89 
90  // Operator
91  operator bool() const {
92  return (numStaticPredictors_ > 0 || numVariablePredictors_ > 0) && vars_.size() > 0;
93  }
94 
95  private:
96  void print(std::ostream &) const override;
97 
98  /// index in the flattened biascoeffs_ for predictor \p jpred and variable \p jvar
99  size_t index(size_t jpred, size_t jvar) const {return jvar*numVariablePredictors_ + jpred;}
100 
101  void initPredictor(const PredictorParametersWrapper &params);
102 
103  /// bias correction coefficients (npredictors x nprimitivevariables)
104  Eigen::VectorXd biascoeffs_;
105 
106  /// bias correction predictors
108  /// predictor names
109  std::vector<std::string> prednames_;
110  /// number of static predictors (i.e. predictors with fixed coefficients all equal to 1)
111  std::size_t numStaticPredictors_;
112  /// number of variable predictors (i.e. predictors with variable coefficients)
114 
115  /// corrected variables names (for now has to be same as "simulated variables")
116  oops::Variables vars_;
117 
118  /// Variables that need to be requested from the model (for computation of predictors)
119  oops::Variables geovars_;
120  /// Diagnostics that need to be requested from the obs operator (for computation of predictors)
121  oops::Variables hdiags_;
122 
123  /// MPI rank, used to determine whether the task should output bias coeffs to a file
124  size_t rank_;
125 };
126 
127 // -----------------------------------------------------------------------------
128 
129 } // namespace ufo
130 
131 #endif // UFO_OBSBIAS_H_
const oops::Variables & correctedVars() const
Return the list of bias-corrected variables.
ObsBiasParameters Parameters_
ObsBias & operator=(const ObsBias &)
Definition: ObsBias.cc:91
size_t rank_
MPI rank, used to determine whether the task should output bias coeffs to a file.
const oops::Variables & requiredVars() const
Eigen::VectorXd biascoeffs_
bias correction coefficients (npredictors x nprimitivevariables)
Predictors predictors_
bias correction predictors
double operator()(size_t jpred, size_t jvar) const
void zero()
Set all variable predictors coefficients to zero (used in the test)
Definition: ObsBias.cc:248
const Eigen::VectorXd & data() const
Return bias correction coefficients (for variable predictors)
oops::Variables geovars_
Variables that need to be requested from the model (for computation of predictors)
ObsBias(ioda::ObsSpace &, const Parameters_ &)
Definition: ObsBias.cc:34
void print(std::ostream &) const override
Definition: ObsBias.cc:261
void initPredictor(const PredictorParametersWrapper &params)
Definition: ObsBias.cc:294
ObsBias & operator+=(const ObsBiasIncrement &)
Definition: ObsBias.cc:84
const Predictors & predictors() const
Return a reference to the vector of all (static and variable) predictors.
const std::vector< std::string > & requiredPredictors() const
std::size_t size() const
static const std::string classname()
std::vector< std::string > prednames_
predictor names
std::size_t numStaticPredictors_
number of static predictors (i.e. predictors with fixed coefficients all equal to 1)
double norm() const
Definition: ObsBias.cc:227
size_t index(size_t jpred, size_t jvar) const
index in the flattened biascoeffs_ for predictor jpred and variable jvar
std::vector< std::shared_ptr< const PredictorBase > > variablePredictors() const
Return the vector of variable predictors.
Definition: ObsBias.cc:254
oops::Variables hdiags_
Diagnostics that need to be requested from the obs operator (for computation of predictors)
void read(const Parameters_ &)
Read bias correction coefficients from file.
Definition: ObsBias.cc:108
oops::Variables vars_
corrected variables names (for now has to be same as "simulated variables")
std::size_t numVariablePredictors_
number of variable predictors (i.e. predictors with variable coefficients)
const oops::Variables & requiredHdiagnostics() const
void write(const Parameters_ &) const
Definition: ObsBias.cc:194
Contains increments to bias correction coefficients.
Parameters influencing the bias correction process.
Contains a polymorphic parameter holding an instance of a subclass of PredictorParametersBase.
Forward declarations.
Definition: ObsAodExt.h:25
Definition: RunCRTM.h:27
std::vector< std::shared_ptr< PredictorBase > > Predictors
Definition: PredictorBase.h:92