8 #ifndef OOPS_ASSIMILATION_GETKFSOLVER_H_
9 #define OOPS_ASSIMILATION_GETKFSOLVER_H_
11 #include <Eigen/Dense>
16 #include "eckit/config/LocalConfiguration.h"
33 #include "oops/util/Logger.h"
34 #include "oops/util/Timer.h"
46 template <
typename MODEL,
typename OBS>
62 static const std::string
classname() {
return "oops::GETKFSolver";}
101 template <
typename MODEL,
typename OBS>
103 const eckit::Configuration & config,
size_t nens)
105 nens_(nens), geometry_(geometry),
106 vertloc_(geometry_, config.getSubConfiguration(
"local ensemble DA.vertical localization")),
107 neig_(vertloc_.neig()), nanal_(neig_*nens_), HZb_(obspaces, nanal_)
113 Log::info() <<
"Multiplicative inflation multCoeff=" <<
114 inflopt.
mult << std::endl;
117 Log::info() <<
"RTPP inflation will be applied with rtppCoeff=" <<
118 inflopt.
rtpp << std::endl;
120 Log::info() <<
"RTPP inflation is not applied rtppCoeff is out of bounds (0,1], rtppCoeff="
121 << inflopt.
rtpp << std::endl;
125 Log::info() <<
"RTPS inflation will be applied with rtpsCoeff=" <<
126 inflopt.
rtps << std::endl;
128 Log::info() <<
"RTPS inflation is not applied rtpsCoeff is out of bounds (0,1], rtpsCoeff="
129 << inflopt.
rtps << std::endl;
138 template <
typename MODEL,
typename OBS>
140 size_t iteration,
bool readFromFile) {
141 util::Timer timer(classname(),
"computeHofX");
148 Log::debug() <<
"Read H(X) from disk" << std::endl;
152 for (
size_t iens = 0; iens < nens_; ++iens) {
153 for (
size_t ieig = 0; ieig < neig_; ++ieig) {
154 ytmp.
read(
"hofxm"+std::to_string(iteration)+
"_"+std::to_string(ieig+1)+
155 "_"+std::to_string(iens+1));
156 HZb_[ii] = ytmp - yb_mean;
157 Log::test() <<
"H(Zx) - ymean for member " << iens+1 <<
" eig "<< ieig+1
158 <<
" :" << std::endl << HZb_[ii] << std::endl;
163 Log::debug() <<
"Computing H(X) online" << std::endl;
169 for (
size_t iens = 0; iens < nens_; ++iens) {
170 vertloc_.modulateIncrement(dx[iens], Ztmp);
171 for (
size_t ieig = 0; ieig < neig_; ++ieig) {
173 tmpState += Ztmp[ieig];
175 HZb_[ii] = tmpObs - yb_mean;
176 tmpObs.
save(
"hofxm"+std::to_string(iteration)+
"_"+std::to_string(ieig+1)+
177 "_"+std::to_string(iens+1));
178 Log::test() <<
"H(Zx) - ymean for member " << iens+1 <<
" eig "<< ieig+1 <<
" :"
179 << std::endl << HZb_[ii] << std::endl;
189 template <
typename MODEL,
typename OBS>
197 util::Timer timer(classname(),
"computeWeights");
199 const int nobsl = dy.
nobs();
204 Eigen::MatrixXf edy_f = edy.cast<
float>();
207 Eigen::MatrixXf eYb_f = eYb.cast<
float>();
209 Eigen::MatrixXd eYb2 = YbOrig.
packEigen();
210 Eigen::MatrixXf eYb2_f = eYb2.cast<
float>();
213 Eigen::MatrixXf eR_f = eR.cast<
float>();
215 Eigen::MatrixXf Wa_f(this->nanal_, this->nens_);
216 Eigen::VectorXf wa_f(this->nanal_);
219 const int getkf_inflation = 0;
223 wa_f.data(), Wa_f.data(),
224 eR_f.data(), nanal_, neig_,
225 getkf_inflation, denkf, getkf);
226 this->Wa_ = Wa_f.cast<
double>();
227 this->wa_ = wa_f.cast<
double>();
232 template <
typename MODEL,
typename OBS>
237 util::Timer timer(classname(),
"applyWeights");
243 std::vector<double> tmp = gptmpl.
getVals();
244 size_t ngp = tmp.size();
245 Eigen::MatrixXd Xb(ngp, nens_);
248 for (
unsigned itime=0; itime < bkg_pert[0].
size(); ++itime) {
251 Xb = vertloc_.modulateIncrement(bkg_pert, i, itime);
254 Eigen::VectorXd xa = Xb*wa_;
255 Eigen::MatrixXd Xa = Xb*Wa_;
259 Xb.resize(ngp, nens_);
260 for (
size_t iens=0; iens < nens_; ++iens) {
262 std::vector<double> tmp1 = gp.
getVals();
263 for (
size_t iv=0; iv < ngp; ++iv) {
264 Xb(iv, iens) = tmp1[iv];
271 Xa = (1-inflopt.
rtpp)*Xa+inflopt.
rtpp*Xb;
275 double eps = DBL_EPSILON;
278 Eigen::ArrayXd asprd = Xa.array().square().rowwise().sum()/(nens_-1);
279 asprd = asprd.sqrt();
280 asprd = (asprd < eps).select(eps, asprd);
283 Eigen::ArrayXd fsprd = Xb.array().square().rowwise().sum()/(nens_-1);
284 fsprd = fsprd.sqrt();
285 fsprd = (fsprd < eps).select(eps, fsprd);
288 Eigen::ArrayXd rtpsInfl = inflopt.
rtps*((fsprd-asprd)/asprd) + 1;
293 Xa.array().colwise() *= rtpsInfl;
297 for (
size_t iens=0; iens < nens_; ++iens) {
298 for (
size_t iv=0; iv < ngp; ++iv) {
299 tmp[iv] = Xa(iv, iens)+xa(iv);
302 ana_pert[iens][itime].setLocal(gptmpl, i);
309 template <
typename MODEL,
typename OBS>
313 util::Timer timer(classname(),
"measurementUpdate");
316 ObsSpaces_ local_obs(this->obspaces_, *i, this->obsconf_);
319 if (local_omb.
nobs() == 0) {
322 this->copyLocalIncrement(bkg_pert, i, ana_pert);
328 ObsErrors_ local_R(this->obsconf_, local_obs);
329 computeWeights(local_omb, local_HZ, local_Yb, local_R);
330 applyWeights(bkg_pert, ana_pert, i);
335 #endif // OOPS_ASSIMILATION_GETKFSOLVER_H_