8 #ifndef OOPS_ASSIMILATION_LOCALENSEMBLESOLVER_H_
9 #define OOPS_ASSIMILATION_LOCALENSEMBLESOLVER_H_
17 #include "eckit/config/Configuration.h"
18 #include "eckit/config/LocalConfiguration.h"
36 #include "oops/util/ConfigFunctions.h"
37 #include "oops/util/Logger.h"
38 #include "oops/util/Timer.h"
43 template <
typename MODEL,
typename OBS>
67 static const std::string
classname() {
return "oops::LocalEnsembleSolver";}
71 const eckit::Configuration & config,
size_t nens);
98 std::unique_ptr<ObsErrors_>
R_;
109 template <
typename MODEL,
typename OBS>
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_)
120 template <
typename MODEL,
typename OBS>
124 const std::vector<util::DateTime> times = xx.
validTimes();
125 const util::Duration flength = times[times.size()-1] - times[0];
130 const util::Duration winlen = obspaces_.windowEnd() - obspaces_.windowStart();
133 std::unique_ptr<PseudoModel_> pseudomodel(
new PseudoModel_(xx, winlen));
134 const Model_ model(std::move(pseudomodel));
136 ModelAux_ moderr(geometry_, eckit::LocalConfiguration());
137 ObsAux_ obsaux(obspaces_, obsconf_);
138 R_.reset(
new ObsErrors_(obsconf_, obspaces_));
144 hofx.
initialize(geometry_, obsaux, *R_, post, config);
145 model.
forecast(init_xx, moderr, flength, post);
151 template <
typename MODEL,
typename OBS>
153 size_t iteration,
bool readFromDisk) {
154 util::Timer timer(classname(),
"computeHofX");
156 ASSERT(ens_xx.
size() == Yb_.size());
158 const size_t nens = ens_xx.
size();
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;
168 R_.reset(
new ObsErrors_(obsconf_, obspaces_));
171 Log::debug() <<
"Computing H(X) online" << std::endl;
172 for (
size_t jj = 0; jj < nens; ++jj) {
173 eckit::LocalConfiguration config;
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));
185 R_->save(
"ObsError");
188 invVarR_.reset(
new Departures_(R_->inverseVariance()));
194 for (
size_t iens = 0; iens < nens; ++iens) {
195 Yb_[iens] = obsens[iens] - yb_mean;
196 Yb_[iens].mask(*invVarR_);
201 omb_ = yobs - yb_mean;
202 omb_.mask(*invVarR_);
210 template <
typename MODEL,
typename OBS>
215 for (
size_t itime=0; itime < bkg_pert[0].
size(); ++itime) {
216 for (
size_t iens=0; iens < bkg_pert.
size(); ++iens) {
218 ana_pert[iens][itime].setLocal(gp,
i);
226 template <
typename MODEL,
typename OBS>
232 const eckit::Configuration &,
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_;
248 template<
class MODEL,
class OBS,
class T>
254 const eckit::Configuration & conf,
size_t nens)
255 {
return new T(obspaces, geometry, conf, nens); }
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.");
268 getMakers()[
name] =
this;
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;
285 jj = getMakers().begin(); jj != getMakers().end(); ++jj) {
286 Log::error() <<
"A " << jj->first <<
" LocalEnsembleSolver" << std::endl;
288 throw std::runtime_error(
id +
" does not exist in local ensemble solver factory.");
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;
Ensemble of Departures (can hold ensemble perturbations in the observation space)
Difference between two observation vectors.
Geometry class used in oops; subclass of interface class interface::Geometry.
Ensemble of 4D increments.
size_t size() const
Accessors.
factory for LETKF solvers
virtual ~LocalEnsembleSolverFactory()=default
Geometry< MODEL > Geometry_
LocalEnsembleSolverFactory(const std::string &)
static std::map< std::string, LocalEnsembleSolverFactory< MODEL, OBS > * > & getMakers()
ObsSpaces< OBS > ObsSpaces_
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.
ObsErrors< OBS > ObsErrors_
const ObsLocalizations_ & obsloc() const
accessor to obs localizations
ObsSpaces< OBS > ObsSpaces_
State4D< MODEL > State4D_
std::unique_ptr< Departures_ > invVarR_
PseudoModelState4D< MODEL > PseudoModel_
Geometry< MODEL > Geometry_
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
Departures< OBS > Departures_
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)
ObsSpaces< OBS > ObsSpaces_
virtual LocalEnsembleSolver< MODEL, OBS > * make(ObsSpaces_ &obspaces, const Geometry_ &geometry, const eckit::Configuration &conf, size_t nens)
Geometry< MODEL > Geometry_
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))
Observations_ mean() const
Compute ensemble mean.
Container for ObsErrors for all observation types that are used in DA.
Container for ObsLocalizations for all observation types that are used in DA.
Computes observation operator (from GeoVaLs), applies bias correction and runs QC filters.
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.
void finalize(Observations_ &)
Computes H(x) from the filled in GeoVaLs.
Control model post processing.
Four dimensional state (vector of 3D States)
const std::vector< util::DateTime > validTimes() const
unsigned int size() const
Accessors.
State class used in oops; subclass of interface class interface::State.
The namespace for the main oops code.