UFO
ObsBiasIncrement.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_OBSBIASINCREMENT_H_
9 #define UFO_OBSBIASINCREMENT_H_
10 
11 #include <string>
12 #include <vector>
13 
14 #include "oops/util/Printable.h"
15 
16 #include "ufo/ObsBiasParameters.h"
18 
19 namespace ioda {
20  class ObsSpace;
21 }
22 
23 namespace ufo {
24  class ObsBias;
25 
26 /// Contains increments to bias correction coefficients
27 class ObsBiasIncrement : public util::Printable {
28  public:
30 
31  ObsBiasIncrement(const ioda::ObsSpace & odb,
32  const Parameters_ & params);
33  ObsBiasIncrement(const ObsBiasIncrement &, const bool = true);
34 
35  // Linear algebra operators
36  void diff(const ObsBias &, const ObsBias &);
37  void zero();
41  ObsBiasIncrement & operator*=(const double);
42  void axpy(const double, const ObsBiasIncrement &);
43  double dot_product_with(const ObsBiasIncrement &) const;
44 
45  // I/O and diagnostics
46  void read(const eckit::Configuration &) {}
47  void write(const eckit::Configuration &) const {}
48  double norm() const;
49 
50  /// Return bias coefficient increments
51  const Eigen::VectorXd & data() const {return biascoeffsinc_;}
52  Eigen::VectorXd & data() {return biascoeffsinc_;}
53 
54  // We could store coefficients in a different order. Then it would
55  // be easier to extract part of the existing vector.
56  /// Return bias coefficient increments for predictor with index \p jpred
57  std::vector<double> coefficients(size_t jpred) const {
58  std::vector<double> coeffs(vars_.size());
59  for (size_t jvar = 0; jvar < vars_.size(); ++jvar) {
60  coeffs[jvar] = biascoeffsinc_(jvar*prednames_.size() + jpred);
61  }
62  return coeffs;
63  }
64  /// Increment bias coeffiecient increments for predictor with index \p jpred
65  /// with \p coeffs
66  void updateCoeff(size_t jpred, const std::vector<double> &coeffs) {
67  for (size_t jj = 0; jj < vars_.size(); ++jj) {
68  biascoeffsinc_[jj*prednames_.size() + jpred] += coeffs[jj];
69  }
70  }
71 
72  // Serialize and deserialize
73  std::size_t serialSize() const {return biascoeffsinc_.size();}
74  void serialize(std::vector<double> &) const;
75  void deserialize(const std::vector<double> &, std::size_t &);
76 
77  // Operator
78  operator bool() const {return biascoeffsinc_.size() > 0;}
79 
80  private:
81  void print(std::ostream &) const;
82 
83  /// Bias coefficient increments
84  Eigen::VectorXd biascoeffsinc_;
85  std::vector<std::string> prednames_;
86 
87  /// Variables that need to be bias-corrected
88  oops::Variables vars_;
89 };
90 
91 // -----------------------------------------------------------------------------
92 
93 } // namespace ufo
94 
95 #endif // UFO_OBSBIASINCREMENT_H_
Contains increments to bias correction coefficients.
void serialize(std::vector< double > &) const
Eigen::VectorXd & data()
ObsBiasIncrement & operator-=(const ObsBiasIncrement &)
double dot_product_with(const ObsBiasIncrement &) const
const Eigen::VectorXd & data() const
Return bias coefficient increments.
void axpy(const double, const ObsBiasIncrement &)
ObsBiasParameters Parameters_
ObsBiasIncrement & operator+=(const ObsBiasIncrement &)
Eigen::VectorXd biascoeffsinc_
Bias coefficient increments.
std::vector< double > coefficients(size_t jpred) const
Return bias coefficient increments for predictor with index jpred.
std::vector< std::string > prednames_
void diff(const ObsBias &, const ObsBias &)
ObsBiasIncrement & operator=(const ObsBiasIncrement &)
std::size_t serialSize() const
void updateCoeff(size_t jpred, const std::vector< double > &coeffs)
void read(const eckit::Configuration &)
oops::Variables vars_
Variables that need to be bias-corrected.
ObsBiasIncrement & operator*=(const double)
ObsBiasIncrement(const ioda::ObsSpace &odb, const Parameters_ &params)
void print(std::ostream &) const
void write(const eckit::Configuration &) const
void deserialize(const std::vector< double > &, std::size_t &)
Parameters influencing the bias correction process.
Forward declarations.
Definition: ObsAodExt.h:25
Definition: RunCRTM.h:27