IODA Bundle
LocalEnsembleSolver.h
Go to the documentation of this file.
1 /*
2  * (C) Copyright 2020 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 OOPS_ASSIMILATION_LOCALENSEMBLESOLVER_H_
9 #define OOPS_ASSIMILATION_LOCALENSEMBLESOLVER_H_
10 
11 #include <map>
12 #include <memory>
13 #include <string>
14 #include <utility>
15 #include <vector>
16 
17 #include "eckit/config/Configuration.h"
18 #include "eckit/config/LocalConfiguration.h"
19 #include "oops/base/Departures.h"
21 #include "oops/base/Geometry.h"
23 #include "oops/base/Model.h"
25 #include "oops/base/ObsEnsemble.h"
26 #include "oops/base/ObsErrors.h"
27 #include "oops/base/Observations.h"
28 #include "oops/base/Observers.h"
30 #include "oops/base/ObsSpaces.h"
31 #include "oops/base/State.h"
36 #include "oops/util/ConfigFunctions.h"
37 #include "oops/util/Logger.h"
38 #include "oops/util/Timer.h"
39 
40 namespace oops {
41 
42 /// \brief Base class for LETKF-type solvers
43 template <typename MODEL, typename OBS>
64  typedef std::vector<std::shared_ptr<ObsData_>> ObsDataVec_;
65 
66  public:
67  static const std::string classname() {return "oops::LocalEnsembleSolver";}
68 
69  /// initialize solver with \p obspaces, \p geometry, full \p config and \p nens ensemble size
70  LocalEnsembleSolver(ObsSpaces_ & obspaces, const Geometry_ & geometry,
71  const eckit::Configuration & config, size_t nens);
72  virtual ~LocalEnsembleSolver() = default;
73 
74  /// computes ensemble H(\p xx), returns mean H(\p xx), saves as hofx \p iteration
75  virtual Observations_ computeHofX(const StateEnsemble4D_ & xx, size_t iteration,
76  bool readFromDisk);
77 
78  /// update background ensemble \p bg to analysis ensemble \p an at a grid point location \p i
79  virtual void measurementUpdate(const IncrementEnsemble4D_ & bg,
80  const GeometryIterator_ & i, IncrementEnsemble4D_ & an) = 0;
81 
82  /// copy \p an[\p i] = \p bg[\p i] (e.g. when there are no local observations to update state)
83  virtual void copyLocalIncrement(const IncrementEnsemble4D_ & bg,
84  const GeometryIterator_ & i, IncrementEnsemble4D_ & an) const;
85 
86  /// compute H(x) based on 4D state \p xx and put the result into \p yy. Also sets up
87  /// R_ based on the QC filters run during H(x)
88  void computeHofX4D(const eckit::Configuration &, const State4D_ &, Observations_ &);
89  /// accessor to obs localizations
90  const ObsLocalizations_ & obsloc() const {return obsloc_;}
91 
92  protected:
93  const Geometry_ & geometry_; ///< Geometry associated with the updated states
94  const ObsSpaces_ & obspaces_; ///< ObsSpaces used in the update
95  Departures_ omb_; ///< obs - mean(H(x)); set in computeHofX method
96  DeparturesEnsemble_ Yb_; ///< ensemble perturbations in the observation space;
97  /// set in computeHofX method
98  std::unique_ptr<ObsErrors_> R_; ///< observation errors, set in computeHofX method
99  std::unique_ptr<Departures_> invVarR_; ///< inverse observation error variance; set in
100  /// computeHofX method
101 
102  private:
103  const eckit::LocalConfiguration obsconf_; // configuration for observations
104  ObsLocalizations_ obsloc_; ///< observation space localization
105 };
106 
107 // -----------------------------------------------------------------------------
108 
109 template <typename MODEL, typename OBS>
111  const Geometry_ & geometry,
112  const eckit::Configuration & config, size_t nens)
113  : geometry_(geometry), obspaces_(obspaces), omb_(obspaces_), Yb_(obspaces_, nens),
114  obsconf_(config, "observations"), obsloc_(obsconf_, obspaces_)
115 {
116 }
117 
118 // -----------------------------------------------------------------------------
119 
120 template <typename MODEL, typename OBS>
121 void LocalEnsembleSolver<MODEL, OBS>::computeHofX4D(const eckit::Configuration & config,
122  const State4D_ & xx, Observations_ & yy) {
123  // compute forecast length from State4D times
124  const std::vector<util::DateTime> times = xx.validTimes();
125  const util::Duration flength = times[times.size()-1] - times[0];
126  // observation window is passed to PseudoModel as the default model time step.
127  // This is required when State4D has a single state: in this case if the default
128  // model time step is not specified, it will be set to zero and no observations
129  // will be processed during the "forecast"
130  const util::Duration winlen = obspaces_.windowEnd() - obspaces_.windowStart();
131 
132  // Setup PseudoModelState4D
133  std::unique_ptr<PseudoModel_> pseudomodel(new PseudoModel_(xx, winlen));
134  const Model_ model(std::move(pseudomodel));
135  // Setup model and obs biases; obs errors
136  ModelAux_ moderr(geometry_, eckit::LocalConfiguration());
137  ObsAux_ obsaux(obspaces_, obsconf_);
138  R_.reset(new ObsErrors_(obsconf_, obspaces_));
139  // Setup and run the model forecast with observers
140  State_ init_xx = xx[0];
142  Observers_ hofx(obspaces_, obsconf_);
143 
144  hofx.initialize(geometry_, obsaux, *R_, post, config);
145  model.forecast(init_xx, moderr, flength, post);
146  hofx.finalize(yy);
147 }
148 
149 // -----------------------------------------------------------------------------
150 
151 template <typename MODEL, typename OBS>
153  size_t iteration, bool readFromDisk) {
154  util::Timer timer(classname(), "computeHofX");
155 
156  ASSERT(ens_xx.size() == Yb_.size());
157 
158  const size_t nens = ens_xx.size();
159  ObsEnsemble_ obsens(obspaces_, nens);
160 
161  if (readFromDisk) {
162  // read hofx from disk
163  Log::debug() << "Read H(X) from disk" << std::endl;
164  for (size_t jj = 0; jj < nens; ++jj) {
165  obsens[jj].read("hofx"+std::to_string(iteration)+"_"+std::to_string(jj+1));
166  Log::test() << "H(x) for member " << jj+1 << ":" << std::endl << obsens[jj] << std::endl;
167  }
168  R_.reset(new ObsErrors_(obsconf_, obspaces_));
169  } else {
170  // compute and save H(x)
171  Log::debug() << "Computing H(X) online" << std::endl;
172  for (size_t jj = 0; jj < nens; ++jj) {
173  eckit::LocalConfiguration config;
174  // never save H(x) (saved explicitly below), save EffectiveQC and EffectiveError
175  // only for the last member
176  config.set("save hofx", false);
177  config.set("save qc", (jj == nens-1));
178  config.set("save obs errors", (jj == nens-1));
179  computeHofX4D(config, ens_xx[jj], obsens[jj]);
180  Log::test() << "H(x) for member " << jj+1 << ":" << std::endl << obsens[jj] << std::endl;
181  obsens[jj].save("hofx"+std::to_string(iteration)+"_"+std::to_string(jj+1));
182  }
183  // QC flags and Obs errors are set to that of the last ensemble member
184  // TODO(someone) combine qc flags from all ensemble members
185  R_->save("ObsError");
186  }
187  // set inverse variances
188  invVarR_.reset(new Departures_(R_->inverseVariance()));
189 
190  // calculate H(x) ensemble mean
191  Observations_ yb_mean(obsens.mean());
192 
193  // calculate H(x) ensemble perturbations
194  for (size_t iens = 0; iens < nens; ++iens) {
195  Yb_[iens] = obsens[iens] - yb_mean;
196  Yb_[iens].mask(*invVarR_);
197  }
198 
199  // calculate obs departures and mask with qc flag
200  Observations_ yobs(obspaces_, "ObsValue");
201  omb_ = yobs - yb_mean;
202  omb_.mask(*invVarR_);
203 
204  // return mean H(x)
205  return yb_mean;
206 }
207 
208 // -----------------------------------------------------------------------------
209 
210 template <typename MODEL, typename OBS>
212  const GeometryIterator_ & i,
213  IncrementEnsemble4D_ & ana_pert) const {
214  // ana_pert[i]=bkg_pert[i]
215  for (size_t itime=0; itime < bkg_pert[0].size(); ++itime) {
216  for (size_t iens=0; iens < bkg_pert.size(); ++iens) {
217  LocalIncrement gp = bkg_pert[iens][itime].getLocal(i);
218  ana_pert[iens][itime].setLocal(gp, i);
219  }
220  }
221 }
222 
223 // =============================================================================
224 
225 /// \brief factory for LETKF solvers
226 template <typename MODEL, typename OBS>
230  public:
231  static std::unique_ptr<LocalEnsembleSolver<MODEL, OBS>> create(ObsSpaces_ &, const Geometry_ &,
232  const eckit::Configuration &,
233  size_t);
234  virtual ~LocalEnsembleSolverFactory() = default;
235  protected:
236  explicit LocalEnsembleSolverFactory(const std::string &);
237  private:
239  const eckit::Configuration &, size_t) = 0;
240  static std::map < std::string, LocalEnsembleSolverFactory<MODEL, OBS> * > & getMakers() {
241  static std::map < std::string, LocalEnsembleSolverFactory<MODEL, OBS> * > makers_;
242  return makers_;
243  }
244 };
245 
246 // -----------------------------------------------------------------------------
247 
248 template<class MODEL, class OBS, class T>
252 
253  virtual LocalEnsembleSolver<MODEL, OBS> * make(ObsSpaces_ & obspaces, const Geometry_ & geometry,
254  const eckit::Configuration & conf, size_t nens)
255  { return new T(obspaces, geometry, conf, nens); }
256  public:
257  explicit LocalEnsembleSolverMaker(const std::string & name)
258  : LocalEnsembleSolverFactory<MODEL, OBS>(name) {}
259 };
260 
261 // =============================================================================
262 
263 template <typename MODEL, typename OBS>
265  if (getMakers().find(name) != getMakers().end()) {
266  throw std::runtime_error(name + " already registered in local ensemble solver factory.");
267  }
268  getMakers()[name] = this;
269 }
270 
271 // -----------------------------------------------------------------------------
272 
273 template <typename MODEL, typename OBS>
274 std::unique_ptr<LocalEnsembleSolver<MODEL, OBS>>
276  const eckit::Configuration & conf, size_t nens) {
277  Log::trace() << "LocalEnsembleSolver<MODEL, OBS>::create starting" << std::endl;
278  const std::string id = conf.getString("local ensemble DA.solver");
279  typename std::map<std::string, LocalEnsembleSolverFactory<MODEL, OBS>*>::iterator
280  jloc = getMakers().find(id);
281  if (jloc == getMakers().end()) {
282  Log::error() << id << " does not exist in local ensemble solver factory." << std::endl;
283  Log::error() << "LETKF solver Factory has " << getMakers().size() << " elements:" << std::endl;
284  for (typename std::map<std::string, LocalEnsembleSolverFactory<MODEL, OBS>*>::const_iterator
285  jj = getMakers().begin(); jj != getMakers().end(); ++jj) {
286  Log::error() << "A " << jj->first << " LocalEnsembleSolver" << std::endl;
287  }
288  throw std::runtime_error(id + " does not exist in local ensemble solver factory.");
289  }
290  std::unique_ptr<LocalEnsembleSolver<MODEL, OBS>>
291  ptr(jloc->second->make(obspaces, geometry, conf, nens));
292  Log::trace() << "LocalEnsembleSolver<MODEL, OBS>::create done" << std::endl;
293  return ptr;
294 }
295 
296 // -----------------------------------------------------------------------------
297 
298 } // namespace oops
299 #endif // OOPS_ASSIMILATION_LOCALENSEMBLESOLVER_H_
program test
Ensemble of Departures (can hold ensemble perturbations in the observation space)
Difference between two observation vectors.
Definition: Departures.h:44
Geometry class used in oops; subclass of interface class interface::Geometry.
Ensemble of 4D increments.
size_t size() const
Accessors.
virtual ~LocalEnsembleSolverFactory()=default
LocalEnsembleSolverFactory(const std::string &)
static std::map< std::string, LocalEnsembleSolverFactory< MODEL, OBS > * > & getMakers()
virtual LocalEnsembleSolver< MODEL, OBS > * make(ObsSpaces_ &, const Geometry_ &, const eckit::Configuration &, size_t)=0
static std::unique_ptr< LocalEnsembleSolver< MODEL, OBS > > create(ObsSpaces_ &, const Geometry_ &, const eckit::Configuration &, size_t)
Base class for LETKF-type solvers.
Observations< OBS > Observations_
GeometryIterator< MODEL > GeometryIterator_
const Geometry_ & geometry_
Geometry associated with the updated states.
const ObsLocalizations_ & obsloc() const
accessor to obs localizations
std::unique_ptr< Departures_ > invVarR_
PseudoModelState4D< MODEL > PseudoModel_
virtual void copyLocalIncrement(const IncrementEnsemble4D_ &bg, const GeometryIterator_ &i, IncrementEnsemble4D_ &an) const
copy an[i] = bg i \a
void computeHofX4D(const eckit::Configuration &, const State4D_ &, Observations_ &)
ObsAuxControls< OBS > ObsAux_
std::unique_ptr< ObsErrors_ > R_
observation errors, set in computeHofX method
ObsLocalizations_ obsloc_
observation space localization
const ObsSpaces_ & obspaces_
ObsSpaces used in the update.
ObsLocalizations< MODEL, OBS > ObsLocalizations_
virtual Observations_ computeHofX(const StateEnsemble4D_ &xx, size_t iteration, bool readFromDisk)
computes ensemble H(xx), returns mean H(xx), saves as hofx iteration
Observers< MODEL, OBS > Observers_
std::vector< std::shared_ptr< ObsData_ > > ObsDataVec_
virtual void measurementUpdate(const IncrementEnsemble4D_ &bg, const GeometryIterator_ &i, IncrementEnsemble4D_ &an)=0
update background ensemble bg to analysis ensemble an at a grid point location i
DeparturesEnsemble< OBS > DeparturesEnsemble_
StateEnsemble4D< MODEL > StateEnsemble4D_
ObsDataVector< OBS, int > ObsData_
Departures_ omb_
obs - mean(H(x)); set in computeHofX method
IncrementEnsemble4D< MODEL > IncrementEnsemble4D_
ObsEnsemble< OBS > ObsEnsemble_
static const std::string classname()
virtual ~LocalEnsembleSolver()=default
const eckit::LocalConfiguration obsconf_
LocalEnsembleSolver(ObsSpaces_ &obspaces, const Geometry_ &geometry, const eckit::Configuration &config, size_t nens)
initialize solver with obspaces, geometry, full config and nens ensemble size
ModelAuxControl< MODEL > ModelAux_
LocalEnsembleSolverMaker(const std::string &name)
virtual LocalEnsembleSolver< MODEL, OBS > * make(ObsSpaces_ &obspaces, const Geometry_ &geometry, const eckit::Configuration &conf, size_t nens)
Abstract nonlinear forecast model used by high level algorithms and applications.
void forecast(State_ &xx, const ModelAux_ &, const util::Duration &len, PostProcessor< State_ > &post) const
Run the forecast from state xx for len time, with post postprocessors Does not need to be implemented...
Ensemble of observations (can hold ensemble of H(x))
Definition: ObsEnsemble.h:21
Observations_ mean() const
Compute ensemble mean.
Definition: ObsEnsemble.h:58
Container for ObsErrors for all observation types that are used in DA.
Definition: ObsErrors.h:32
Container for ObsLocalizations for all observation types that are used in DA.
Observations Class.
Definition: Observations.h:35
Computes observation operator (from GeoVaLs), applies bias correction and runs QC filters.
Definition: Observers.h:36
void initialize(const Geometry_ &, const ObsAuxCtrls_ &, ObsErrors_ &, PostProc_ &, const eckit::Configuration &=eckit::LocalConfiguration())
Initializes variables, obs bias, obs filters (could be different for different iterations.
Definition: Observers.h:98
void finalize(Observations_ &)
Computes H(x) from the filled in GeoVaLs.
Definition: Observers.h:115
Control model post processing.
Definition: PostProcessor.h:30
Four dimensional state (vector of 3D States)
Definition: State4D.h:29
const std::vector< util::DateTime > validTimes() const
Definition: State4D.h:107
Ensemble of 4D states.
unsigned int size() const
Accessors.
State class used in oops; subclass of interface class interface::State.
The namespace for the main oops code.