OOPS
ObsBiasIncrement.cc
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 #include "model/ObsBiasIncrement.h"
12 
13 #include <cmath>
14 #include <iostream>
15 #include <sstream>
16 #include <string>
17 
18 #include "eckit/config/Configuration.h"
19 #include "model/ObsBias.h"
21 #include "oops/util/Logger.h"
22 
23 // -----------------------------------------------------------------------------
24 namespace qg {
25 // -----------------------------------------------------------------------------
26 ObsBiasIncrement::ObsBiasIncrement(const ObsSpaceQG &, const eckit::Configuration & conf)
27  : bias_(ObsBias::ntypes, 0.0), active_(ObsBias::ntypes, false)
28 {
29  if (conf.has("obs bias error")) {
30  const eckit::LocalConfiguration covconf(conf, "obs bias error");
31  active_[0] = covconf.has("stream");
32  active_[1] = covconf.has("uwind");
33  active_[2] = covconf.has("vwind");
34  active_[3] = covconf.has("wspeed");
35  }
36  bool on = false;
37  std::string strn = "";
38  for (unsigned int jj = 0; jj < ObsBias::ntypes; ++jj) {
39  if (jj > 0) strn += ", ";
40  if (active_[jj]) {
41  strn += "on";
42  on = true;
43  } else {
44  strn += "off";
45  }
46  }
47  if (on) {oops::Log::trace() << "ObsBiasIncrement created : " << strn << std::endl;}
48 }
49 // -----------------------------------------------------------------------------
51  const bool copy)
52  : bias_(ObsBias::ntypes, 0.0), active_(other.active_)
53 {
54  if (copy) {
55  for (unsigned int jj = 0; jj < ObsBias::ntypes; ++jj) bias_[jj] = other.bias_[jj];
56  }
57  this->makePassive();
58 }
59 // -----------------------------------------------------------------------------
61  const eckit::Configuration &)
62  : bias_(ObsBias::ntypes, 0.0), active_(other.active_)
63 {
64  for (unsigned int jj = 0; jj < ObsBias::ntypes; ++jj) bias_[jj] = other.bias_[jj];
65  this->makePassive();
66 }
67 // -----------------------------------------------------------------------------
69  for (unsigned int jj = 0; jj < ObsBias::ntypes; ++jj) {
70  if (!active_[jj]) bias_[jj] = 0.0;
71  }
72 }
73 // -----------------------------------------------------------------------------
74 void ObsBiasIncrement::diff(const ObsBias & b1, const ObsBias & b2) {
75  for (unsigned int jj = 0; jj < ObsBias::ntypes; ++jj) {
76  bias_[jj] = b1[jj] - b2[jj];
77  }
78  this->makePassive();
79 }
80 // -----------------------------------------------------------------------------
82  for (unsigned int jj = 0; jj < ObsBias::ntypes; ++jj) bias_[jj] = 0.0;
83 }
84 // -----------------------------------------------------------------------------
86  for (unsigned int jj = 0; jj < ObsBias::ntypes; ++jj) bias_[jj] = rhs.bias_[jj];
87  this->makePassive();
88  return *this;
89 }
90 // -----------------------------------------------------------------------------
92  for (unsigned int jj = 0; jj < ObsBias::ntypes; ++jj) bias_[jj] += rhs.bias_[jj];
93  this->makePassive();
94  return *this;
95 }
96 // -----------------------------------------------------------------------------
98  for (unsigned int jj = 0; jj < ObsBias::ntypes; ++jj) bias_[jj] -= rhs.bias_[jj];
99  this->makePassive();
100  return *this;
101 }
102 // -----------------------------------------------------------------------------
104  for (unsigned int jj = 0; jj < ObsBias::ntypes; ++jj) bias_[jj] *= fact;
105  this->makePassive();
106  return *this;
107 }
108 // -----------------------------------------------------------------------------
109 void ObsBiasIncrement::axpy(const double fact, const ObsBiasIncrement & rhs) {
110  for (unsigned int jj = 0; jj < ObsBias::ntypes; ++jj) bias_[jj] += fact * rhs.bias_[jj];
111  this->makePassive();
112 }
113 // -----------------------------------------------------------------------------
115  double zz = 0.0;
116  for (unsigned int jj = 0; jj < ObsBias::ntypes; ++jj) {
117  if (active_[jj]) zz += bias_[jj] * rhs.bias_[jj];
118  }
119  return zz;
120 }
121 // -----------------------------------------------------------------------------
122 double ObsBiasIncrement::norm() const {
123  double zz = 0.0;
124  int ii = 0;
125  for (unsigned int jj = 0; jj < ObsBias::ntypes; ++jj) {
126  if (active_[jj]) {
127  zz += bias_[jj] * bias_[jj];
128  ++ii;
129  }
130  }
131  if (ii > 0) zz = std::sqrt(zz/ii);
132  return zz;
133 }
134 // -----------------------------------------------------------------------------
136  size_t nn = bias_.size();
137  return nn;
138 }
139 // -----------------------------------------------------------------------------
140 void ObsBiasIncrement::serialize(std::vector<double> & vect) const {
141  vect.insert(vect.end(), bias_.begin(), bias_.end());
142  oops::Log::trace() << "ObsBiasIncrement::serialize done" << std::endl;
143 }
144 // -----------------------------------------------------------------------------
145 void ObsBiasIncrement::deserialize(const std::vector<double> & vect, size_t & index) {
146  for (unsigned int jj = 0; jj < bias_.size(); ++jj) {
147  bias_[jj] = vect[index];
148  ++index;
149  }
150  oops::Log::trace() << "ObsBiasIncrement::deserialize done" << std::endl;
151 }
152 // -----------------------------------------------------------------------------
153 void ObsBiasIncrement::print(std::ostream & os) const {
154  bool on = false;
155  std::string strn = "";
156  for (unsigned int jj = 0; jj < ObsBias::ntypes; ++jj) {
157  if (jj > 0) strn += ", ";
158  if (active_[jj]) {
159  on = true;
160  std::ostringstream strs;
161  strs << bias_[jj];
162  strn += strs.str();
163  } else {
164  strn += "0.0";
165  }
166  }
167  if (on) os << std::endl << "ObsBiasIncrement = " << strn;
168 }
169 // -----------------------------------------------------------------------------
170 } // namespace qg
qg::ObsBiasIncrement::dot_product_with
double dot_product_with(const ObsBiasIncrement &) const
Definition: ObsBiasIncrement.cc:114
qg
The namespace for the qg model.
Definition: qg/model/AnalyticInit.cc:13
ObsBiasIncrement.h
qg::ObsBiasIncrement::makePassive
void makePassive()
Definition: ObsBiasIncrement.cc:68
qg::ObsBiasIncrement::ObsBiasIncrement
ObsBiasIncrement()
Constructor, destructor.
qg::ObsBiasIncrement::operator-=
ObsBiasIncrement & operator-=(const ObsBiasIncrement &)
Definition: ObsBiasIncrement.cc:97
qg::ObsBias
Class to handle observation bias parameters.
Definition: qg/model/ObsBias.h:37
qg::ObsBiasIncrement::zero
void zero()
Definition: ObsBiasIncrement.cc:81
qg::ObsBiasIncrement::active_
std::vector< bool > active_
Definition: ObsBiasIncrement.h:75
qg::ObsSpaceQG
ObsSpace for QG model.
Definition: ObsSpaceQG.h:44
qg::ObsBiasIncrement::serialSize
size_t serialSize() const override
Serialization.
Definition: ObsBiasIncrement.cc:135
ObsBias.h
qg::ObsBiasIncrement::print
void print(std::ostream &) const override
Definition: ObsBiasIncrement.cc:153
qg::ObsBiasIncrement::operator=
ObsBiasIncrement & operator=(const ObsBiasIncrement &)
Definition: ObsBiasIncrement.cc:85
ObsBiasCovariance.h
qg::ObsBias::ntypes
static const unsigned int ntypes
Definition: qg/model/ObsBias.h:39
qg::ObsBiasIncrement::bias_
std::vector< double > bias_
Definition: ObsBiasIncrement.h:74
qg::ObsBiasIncrement
Definition: ObsBiasIncrement.h:31
qg::ObsBiasIncrement::diff
void diff(const ObsBias &, const ObsBias &)
Linear algebra operators.
Definition: ObsBiasIncrement.cc:74
qg::ObsBiasIncrement::serialize
void serialize(std::vector< double > &) const override
Definition: ObsBiasIncrement.cc:140
qg::ObsBiasIncrement::operator+=
ObsBiasIncrement & operator+=(const ObsBiasIncrement &)
Definition: ObsBiasIncrement.cc:91
qg::ObsBiasIncrement::deserialize
void deserialize(const std::vector< double > &, size_t &) override
Definition: ObsBiasIncrement.cc:145
qg::ObsBiasIncrement::axpy
void axpy(const double, const ObsBiasIncrement &)
Definition: ObsBiasIncrement.cc:109
qg::ObsBiasIncrement::norm
double norm() const
Definition: ObsBiasIncrement.cc:122
qg::ObsBiasIncrement::operator*=
ObsBiasIncrement & operator*=(const double)
Definition: ObsBiasIncrement.cc:103