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>
59 static const std::string
classname() {
return "oops::LETKFSolver";}
73 virtual void computeWeights(
const Eigen::VectorXd & omb,
const Eigen::MatrixXd & Yb,
74 const Eigen::VectorXd & invvarR);
94 template <
typename MODEL,
typename OBS>
96 const eckit::Configuration & config,
size_t nens,
104 Log::info() <<
"Using EIGEN implementation of LETKF" << std::endl;
106 Log::info() <<
"Multiplicative inflation multCoeff=" <<
107 inflopt.
mult << std::endl;
110 Log::info() <<
"RTPP inflation will be applied with rtppCoeff=" <<
111 inflopt.
rtpp << std::endl;
113 Log::info() <<
"RTPP inflation is not applied rtppCoeff is out of bounds (0,1], rtppCoeff="
114 << inflopt.
rtpp << std::endl;
118 Log::info() <<
"RTPS inflation will be applied with rtpsCoeff=" <<
119 inflopt.
rtps << std::endl;
121 Log::info() <<
"RTPS inflation is not applied rtpsCoeff is out of bounds (0,1], rtpsCoeff="
122 << inflopt.
rtps << std::endl;
136 template <
typename MODEL,
typename OBS>
140 util::Timer timer(classname(),
"measurementUpdate");
145 this->obsloc().computeLocalization(i, locvector);
146 locvector.
mask(*(this->invVarR_));
147 Eigen::VectorXd local_omb_vec = this->omb_.packEigen(locvector);
149 if (local_omb_vec.size() == 0) {
152 this->copyLocalIncrement(bkg_pert, i, ana_pert);
156 Eigen::MatrixXd local_Yb_mat = this->Yb_.packEigen(locvector);
158 Eigen::VectorXd local_invVarR_vec = this->invVarR_->packEigen(locvector);
160 Eigen::VectorXd localization = locvector.
packEigen(locvector);
161 local_invVarR_vec.array() *= localization.array();
162 computeWeights(local_omb_vec, local_Yb_mat, local_invVarR_vec);
163 applyWeights(bkg_pert, ana_pert, i);
169 template <
typename MODEL,
typename OBS>
171 const Eigen::MatrixXd & Yb,
172 const Eigen::VectorXd & diagInvR ) {
176 util::Timer timer(classname(),
"computeWeights");
182 double infl = inflopt.
mult;
183 Eigen::MatrixXd work = Yb*(diagInvR.asDiagonal()*Yb.transpose());
184 work.diagonal() += Eigen::VectorXd::Constant(nens_, (nens_-1)/infl);
187 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> es(work);
188 eival_ = es.eigenvalues().real();
189 eivec_ = es.eigenvectors().real();
192 work = eivec_ * eival_.cwiseInverse().asDiagonal() * eivec_.transpose();
196 * ((nens_-1) * eival_.array().inverse()).sqrt().matrix().asDiagonal()
197 * eivec_.transpose();
200 wa_ = work * (Yb * (diagInvR.asDiagonal()*dy));
205 template <
typename MODEL,
typename OBS>
210 util::Timer timer(classname(),
"applyWeights");
215 std::vector<double> tmp1 = gptmpl.
getVals();
216 size_t ngp = tmp1.size();
219 for (
size_t itime=0; itime < bkg_pert[0].
size(); ++itime) {
221 Eigen::MatrixXd Xb(ngp, nens_);
223 for (
size_t iens=0; iens < nens_; ++iens) {
225 std::vector<double> tmp = gp.
getVals();
226 for (
size_t iv=0; iv < ngp; ++iv) {
227 Xb(iv, iens) = tmp[iv];
232 Eigen::VectorXd xa = Xb*wa_;
233 Eigen::MatrixXd Xa = Xb*Wa_;
237 Xa = (1-inflopt.
rtpp)*Xa+inflopt.
rtpp*Xb;
241 double eps = DBL_EPSILON;
244 Eigen::ArrayXd asprd = Xa.array().square().rowwise().sum()/(nens_-1);
245 asprd = asprd.sqrt();
246 asprd = (asprd < eps).select(eps, asprd);
249 Eigen::ArrayXd fsprd = Xb.array().square().rowwise().sum()/(nens_-1);
250 fsprd = fsprd.sqrt();
251 fsprd = (fsprd < eps).select(eps, fsprd);
254 Eigen::ArrayXd rtpsInfl = inflopt.
rtps*((fsprd-asprd)/asprd) + 1;
259 Xa.array().colwise() *= rtpsInfl;
264 for (
size_t iens=0; iens < nens_; ++iens) {
265 for (
size_t iv=0; iv < ngp; ++iv) {
266 tmp1[iv] = Xa(iv, iens)+xa(iv);
269 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_
LETKFSolver(ObsSpaces_ &, const Geometry_ &, const eckit::Configuration &, size_t, const State4D_ &)
State4D< MODEL > State4D_
LETKFSolverParameters options_
Geometry< MODEL > Geometry_
IncrementEnsemble4D< MODEL > IncrementEnsemble4D_
DeparturesEnsemble< OBS > DeparturesEnsemble_
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.
Four dimensional state (vector of 3D States)
The namespace for the main oops code.