9 #ifndef OOPS_ASSIMILATION_LETKFSOLVER_H_
10 #define OOPS_ASSIMILATION_LETKFSOLVER_H_
12 #include <Eigen/Dense>
18 #include "eckit/config/LocalConfiguration.h"
30 #include "oops/util/Logger.h"
31 #include "oops/util/Timer.h"
46 template <
typename MODEL,
typename OBS>
58 static const std::string
classname() {
return "oops::LETKFSolver";}
71 virtual void computeWeights(
const Eigen::VectorXd & omb,
const Eigen::MatrixXd & Yb,
72 const Eigen::VectorXd & invvarR);
92 template <
typename MODEL,
typename OBS>
94 const eckit::Configuration & config,
size_t nens)
100 Log::info() <<
"Using EIGEN implementation of LETKF" << std::endl;
102 Log::info() <<
"Multiplicative inflation multCoeff=" <<
103 inflopt.
mult << std::endl;
106 Log::info() <<
"RTPP inflation will be applied with rtppCoeff=" <<
107 inflopt.
rtpp << std::endl;
109 Log::info() <<
"RTPP inflation is not applied rtppCoeff is out of bounds (0,1], rtppCoeff="
110 << inflopt.
rtpp << std::endl;
114 Log::info() <<
"RTPS inflation will be applied with rtpsCoeff=" <<
115 inflopt.
rtps << std::endl;
117 Log::info() <<
"RTPS inflation is not applied rtpsCoeff is out of bounds (0,1], rtpsCoeff="
118 << inflopt.
rtps << std::endl;
132 template <
typename MODEL,
typename OBS>
136 util::Timer timer(classname(),
"measurementUpdate");
141 this->obsloc().computeLocalization(
i, locvector);
142 locvector.
mask(*(this->invVarR_));
143 Eigen::VectorXd local_omb_vec = this->omb_.packEigen(locvector);
145 if (local_omb_vec.size() == 0) {
148 this->copyLocalIncrement(bkg_pert,
i, ana_pert);
152 Eigen::MatrixXd local_Yb_mat = this->Yb_.packEigen(locvector);
154 Eigen::VectorXd local_invVarR_vec = this->invVarR_->packEigen(locvector);
156 Eigen::VectorXd localization = locvector.
packEigen(locvector);
157 local_invVarR_vec.array() *= localization.array();
158 computeWeights(local_omb_vec, local_Yb_mat, local_invVarR_vec);
159 applyWeights(bkg_pert, ana_pert,
i);
165 template <
typename MODEL,
typename OBS>
167 const Eigen::MatrixXd & Yb,
168 const Eigen::VectorXd & diagInvR ) {
172 util::Timer timer(classname(),
"computeWeights");
178 double infl = inflopt.
mult;
179 Eigen::MatrixXd work = Yb*(diagInvR.asDiagonal()*Yb.transpose());
180 work.diagonal() += Eigen::VectorXd::Constant(nens_, (nens_-1)/infl);
183 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> es(work);
184 eival_ = es.eigenvalues().real();
185 eivec_ = es.eigenvectors().real();
188 work = eivec_ * eival_.cwiseInverse().asDiagonal() * eivec_.transpose();
192 * ((nens_-1) * eival_.array().inverse()).sqrt().matrix().asDiagonal()
193 * eivec_.transpose();
196 wa_ = work * (Yb * (diagInvR.asDiagonal()*dy));
201 template <
typename MODEL,
typename OBS>
206 util::Timer timer(classname(),
"applyWeights");
211 std::vector<double> tmp1 = gptmpl.
getVals();
212 size_t ngp = tmp1.size();
215 for (
size_t itime=0; itime < bkg_pert[0].
size(); ++itime) {
217 Eigen::MatrixXd Xb(ngp, nens_);
219 for (
size_t iens=0; iens < nens_; ++iens) {
221 std::vector<double> tmp = gp.
getVals();
222 for (
size_t iv=0; iv < ngp; ++iv) {
223 Xb(iv, iens) = tmp[iv];
228 Eigen::VectorXd xa = Xb*wa_;
229 Eigen::MatrixXd Xa = Xb*Wa_;
233 Xa = (1-inflopt.
rtpp)*Xa+inflopt.
rtpp*Xb;
237 double eps = DBL_EPSILON;
240 Eigen::ArrayXd asprd = Xa.array().square().rowwise().sum()/(nens_-1);
241 asprd = asprd.sqrt();
242 asprd = (asprd < eps).select(eps, asprd);
245 Eigen::ArrayXd fsprd = Xb.array().square().rowwise().sum()/(nens_-1);
246 fsprd = fsprd.sqrt();
247 fsprd = (fsprd < eps).select(eps, fsprd);
250 Eigen::ArrayXd rtpsInfl = inflopt.
rtps*((fsprd-asprd)/asprd) + 1;
255 Xa.array().colwise() *= rtpsInfl;
260 for (
size_t iens=0; iens < nens_; ++iens) {
261 for (
size_t iv=0; iv < ngp; ++iv) {
262 tmp1[iv] = Xa(iv, iens)+xa(iv);
265 ana_pert[iens][itime].setLocal(gptmpl,
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)
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
Local Ensemble Tranform Kalman Filter solver.
Departures< OBS > Departures_
void measurementUpdate(const IncrementEnsemble4D_ &, const GeometryIterator_ &, IncrementEnsemble4D_ &) override
KF update + posterior inflation at a grid point location (GeometryIterator_)
ObsLocalizations< MODEL, OBS > ObsLocalizations_
LETKFSolverParameters options_
Geometry< MODEL > Geometry_
IncrementEnsemble4D< MODEL > IncrementEnsemble4D_
DeparturesEnsemble< OBS > DeparturesEnsemble_
LETKFSolver(ObsSpaces_ &, const Geometry_ &, const eckit::Configuration &, size_t)
ObsErrors< OBS > ObsErrors_
virtual void applyWeights(const IncrementEnsemble4D_ &, IncrementEnsemble4D_ &, const GeometryIterator_ &)
Applies weights and adds posterior inflation.
virtual void computeWeights(const Eigen::VectorXd &omb, const Eigen::MatrixXd &Yb, const Eigen::VectorXd &invvarR)
GeometryIterator< MODEL > GeometryIterator_
ObsSpaces< OBS > ObsSpaces_
static const std::string classname()
Parameter< LETKFInflationParameters > infl
Base class for LETKF-type solvers.
void setVals(std::vector< double > &)
const std::vector< double > & getVals() const
Container for ObsErrors for all observation types that are used in DA.
Container for ObsLocalizations for all observation types that are used in DA.
The namespace for the main oops code.