8 #ifndef OOPS_ASSIMILATION_GETKFSOLVER_H_
9 #define OOPS_ASSIMILATION_GETKFSOLVER_H_
11 #include <Eigen/Dense>
18 #include "eckit/config/LocalConfiguration.h"
34 #include "oops/util/ConfigFunctions.h"
35 #include "oops/util/Logger.h"
36 #include "oops/util/Timer.h"
48 template <
typename MODEL,
typename OBS>
63 static const std::string
classname() {
return "oops::GETKFSolver";}
83 void computeWeights(
const Eigen::VectorXd & omb,
const Eigen::MatrixXd & Yb,
84 const Eigen::MatrixXd & YbOrig,
const Eigen::VectorXd & invvarR);
107 template <
typename MODEL,
typename OBS>
109 const eckit::Configuration & config,
size_t nens,
112 nens_(nens), geometry_(geometry),
113 vertloc_(config.getSubConfiguration(
"local ensemble DA.vertical localization"), xbmean[0]),
114 neig_(vertloc_.neig()), nanal_(neig_*nens_), HZb_(obspaces, nanal_)
120 Log::info() <<
"Multiplicative inflation multCoeff=" <<
121 inflopt.
mult << std::endl;
124 Log::info() <<
"RTPP inflation will be applied with rtppCoeff=" <<
125 inflopt.
rtpp << std::endl;
127 Log::info() <<
"RTPP inflation is not applied rtppCoeff is out of bounds (0,1], rtppCoeff="
128 << inflopt.
rtpp << std::endl;
132 Log::info() <<
"RTPS inflation will be applied with rtpsCoeff=" <<
133 inflopt.
rtps << std::endl;
135 Log::info() <<
"RTPS inflation is not applied rtpsCoeff is out of bounds (0,1], rtpsCoeff="
136 << inflopt.
rtps << std::endl;
145 template <
typename MODEL,
typename OBS>
147 size_t iteration,
bool readFromFile) {
148 util::Timer timer(classname(),
"computeHofX");
155 Log::debug() <<
"Read H(X) from disk" << std::endl;
159 for (
size_t iens = 0; iens < nens_; ++iens) {
160 for (
size_t ieig = 0; ieig < neig_; ++ieig) {
161 ytmp.
read(
"hofxm"+std::to_string(iteration)+
"_"+std::to_string(ieig+1)+
162 "_"+std::to_string(iens+1));
163 HZb_[ii] = ytmp - yb_mean;
164 Log::test() <<
"H(Zx) - ymean for member " << iens+1 <<
" eig "<< ieig+1
165 <<
" :" << std::endl << HZb_[ii] << std::endl;
170 Log::debug() <<
"Computing H(X) online" << std::endl;
176 for (
size_t iens = 0; iens < nens_; ++iens) {
177 vertloc_.modulateIncrement(dx[iens], Ztmp);
178 for (
size_t ieig = 0; ieig < neig_; ++ieig) {
180 tmpState += Ztmp[ieig];
182 eckit::LocalConfiguration config;
183 config.set(
"save hofx",
false);
184 config.set(
"save qc",
false);
185 config.set(
"save obs errors",
false);
186 this->computeHofX4D(config, tmpState, tmpObs);
187 HZb_[ii] = tmpObs - yb_mean;
188 tmpObs.
save(
"hofxm"+std::to_string(iteration)+
"_"+std::to_string(ieig+1)+
189 "_"+std::to_string(iens+1));
190 Log::test() <<
"H(Zx) - ymean for member " << iens+1 <<
" eig "<< ieig+1 <<
" :"
191 << std::endl << HZb_[ii] << std::endl;
201 template <
typename MODEL,
typename OBS>
203 const Eigen::MatrixXd & Yb,
204 const Eigen::MatrixXd & YbOrig,
205 const Eigen::VectorXd & R_invvar) {
209 util::Timer timer(classname(),
"computeWeights");
211 const int nobsl = dy.size();
214 Eigen::VectorXf dy_f = dy.cast<
float>();
215 Eigen::MatrixXf Yb_f = Yb.cast<
float>();
216 Eigen::MatrixXf YbOrig_f = YbOrig.cast<
float>();
217 Eigen::VectorXf R_invvar_f = R_invvar.cast<
float>();
219 Eigen::MatrixXf Wa_f(nanal_, this->nens_);
220 Eigen::VectorXf wa_f(nanal_);
223 const int getkf_inflation = 0;
227 wa_f.data(), Wa_f.data(),
228 R_invvar_f.data(), nanal_, neig_,
229 getkf_inflation, denkf, getkf);
230 this->Wa_ = Wa_f.cast<
double>();
231 this->wa_ = wa_f.cast<
double>();
236 template <
typename MODEL,
typename OBS>
241 util::Timer timer(classname(),
"applyWeights");
247 std::vector<double> tmp = gptmpl.
getVals();
248 size_t ngp = tmp.size();
249 Eigen::MatrixXd Xb(ngp, nens_);
252 for (
unsigned itime=0; itime < bkg_pert[0].
size(); ++itime) {
255 Xb = vertloc_.modulateIncrement(bkg_pert, i, itime);
258 Eigen::VectorXd xa = Xb*wa_;
259 Eigen::MatrixXd Xa = Xb*Wa_;
263 Xb.resize(ngp, nens_);
264 for (
size_t iens=0; iens < nens_; ++iens) {
266 std::vector<double> tmp1 = gp.
getVals();
267 for (
size_t iv=0; iv < ngp; ++iv) {
268 Xb(iv, iens) = tmp1[iv];
275 Xa = (1-inflopt.
rtpp)*Xa+inflopt.
rtpp*Xb;
279 double eps = DBL_EPSILON;
282 Eigen::ArrayXd asprd = Xa.array().square().rowwise().sum()/(nens_-1);
283 asprd = asprd.sqrt();
284 asprd = (asprd < eps).select(eps, asprd);
287 Eigen::ArrayXd fsprd = Xb.array().square().rowwise().sum()/(nens_-1);
288 fsprd = fsprd.sqrt();
289 fsprd = (fsprd < eps).select(eps, fsprd);
292 Eigen::ArrayXd rtpsInfl = inflopt.
rtps*((fsprd-asprd)/asprd) + 1;
297 Xa.array().colwise() *= rtpsInfl;
301 for (
size_t iens=0; iens < nens_; ++iens) {
302 for (
size_t iv=0; iv < ngp; ++iv) {
303 tmp[iv] = Xa(iv, iens)+xa(iv);
306 ana_pert[iens][itime].setLocal(gptmpl, i);
313 template <
typename MODEL,
typename OBS>
317 util::Timer timer(classname(),
"measurementUpdate");
322 this->obsloc().computeLocalization(i, locvector);
323 locvector.
mask(*(this->invVarR_));
324 Eigen::VectorXd local_omb_vec = this->omb_.packEigen(locvector);
326 if (local_omb_vec.size() == 0) {
329 this->copyLocalIncrement(bkg_pert, i, ana_pert);
333 Eigen::MatrixXd local_Yb_mat = this->Yb_.packEigen(locvector);
334 Eigen::MatrixXd local_HZ_mat = this->HZb_.packEigen(locvector);
336 Eigen::VectorXd local_invVarR_vec = this->invVarR_->packEigen(locvector);
338 Eigen::VectorXd localization = locvector.
packEigen(locvector);
339 local_invVarR_vec.array() *= localization.array();
340 computeWeights(local_omb_vec, local_HZ_mat, local_Yb_mat, local_invVarR_vec);
341 applyWeights(bkg_pert, ana_pert, i);
Ensemble of Departures (can hold ensemble perturbations in the observation space)
Difference between two observation vectors.
void mask(ObsDataVec_< int >)
Mask out departures where the passed in qc flags are > 0.
Eigen::VectorXd packEigen(const Departures &) const
Pack departures in an Eigen vector (excluding departures that are masked out)
An implementation of the GETKF from Lei 2018 JAMES.
Departures< OBS > Departures_
void applyWeights(const IncrementEnsemble4D_ &, IncrementEnsemble4D_ &, const GeometryIterator_ &)
Applies weights and adds posterior inflation.
VerticalLocEV< MODEL > VerticalLocEV_
GeometryIterator< MODEL > GeometryIterator_
GETKFSolver(ObsSpaces_ &, const Geometry_ &, const eckit::Configuration &, size_t, const State4D_ &)
static const std::string classname()
DeparturesEnsemble< OBS > DeparturesEnsemble_
const Geometry_ & geometry_
IncrementEnsemble4D< MODEL > IncrementEnsemble4D_
ObsEnsemble< OBS > ObsEnsemble_
Geometry< MODEL > Geometry_
void computeWeights(const Eigen::VectorXd &omb, const Eigen::MatrixXd &Yb, const Eigen::MatrixXd &YbOrig, const Eigen::VectorXd &invvarR)
StateEnsemble4D< MODEL > StateEnsemble4D_
Observations_ computeHofX(const StateEnsemble4D_ &, size_t, bool) override
computes ensemble H(xx), returns mean H(xx), saves as hofx iteration
void measurementUpdate(const IncrementEnsemble4D_ &, const GeometryIterator_ &, IncrementEnsemble4D_ &) override
entire KF update (computeWeights+applyWeights) for a grid point GeometryIterator_
State4D< MODEL > State4D_
ObsSpaces< OBS > ObsSpaces_
Observations< OBS > Observations_
LETKFSolverParameters options_
Geometry class used in oops; subclass of interface class interface::Geometry.
Ensemble of 4D increments.
size_t size() const
Accessors.
Parameters for LETKF inflation.
double rtpsInflMax() const
double rtpsInflMin() const
Parameter< LETKFInflationParameters > infl
Base class for LETKF-type solvers.
virtual Observations_ computeHofX(const StateEnsemble4D_ &xx, size_t iteration, bool readFromDisk)
computes ensemble H(xx), returns mean H(xx), saves as hofx iteration
void setVals(std::vector< double > &)
const std::vector< double > & getVals() const
Ensemble of observations (can hold ensemble of H(x))
void save(const std::string &) const
Save/read observations values.
void read(const std::string &)
Four dimensional state (vector of 3D States)
State4D_ mean() const
calculate ensemble mean
The namespace for the main oops code.
void letkf_core_f90(const int &, const float *, const float *, const float *, float *, float *, const float *, const int &, const int &, const int &, const int &, const int &)