IODA Bundle
|
Local Ensemble Tranform Kalman Filter solver. More...
#include <LETKFSolver.h>
Public Member Functions | |
LETKFSolver (ObsSpaces_ &, const Geometry_ &, const eckit::Configuration &, size_t) | |
void | measurementUpdate (const IncrementEnsemble4D_ &, const GeometryIterator_ &, IncrementEnsemble4D_ &) override |
KF update + posterior inflation at a grid point location (GeometryIterator_) More... | |
![]() | |
LocalEnsembleSolver (ObsSpaces_ &obspaces, const Geometry_ &geometry, const eckit::Configuration &config, size_t nens) | |
initialize solver with obspaces , geometry , full config and nens ensemble size More... | |
virtual | ~LocalEnsembleSolver ()=default |
virtual Observations_ | computeHofX (const StateEnsemble4D_ &xx, size_t iteration, bool readFromDisk) |
computes ensemble H(xx ), returns mean H(xx ), saves as hofx iteration More... | |
virtual void | copyLocalIncrement (const IncrementEnsemble4D_ &bg, const GeometryIterator_ &i, IncrementEnsemble4D_ &an) const |
copy an [i ] = bg i \a More... | |
void | computeHofX4D (const eckit::Configuration &, const State4D_ &, Observations_ &) |
const ObsLocalizations_ & | obsloc () const |
accessor to obs localizations More... | |
Static Public Member Functions | |
static const std::string | classname () |
![]() | |
static const std::string | classname () |
Protected Member Functions | |
virtual void | computeWeights (const Eigen::VectorXd &omb, const Eigen::MatrixXd &Yb, const Eigen::VectorXd &invvarR) |
virtual void | applyWeights (const IncrementEnsemble4D_ &, IncrementEnsemble4D_ &, const GeometryIterator_ &) |
Applies weights and adds posterior inflation. More... | |
Protected Attributes | |
LETKFSolverParameters | options_ |
Eigen::MatrixXd | Wa_ |
Eigen::VectorXd | wa_ |
Eigen::VectorXd | eival_ |
Eigen::MatrixXd | eivec_ |
const size_t | nens_ |
![]() | |
const Geometry_ & | geometry_ |
Geometry associated with the updated states. More... | |
const ObsSpaces_ & | obspaces_ |
ObsSpaces used in the update. More... | |
Departures_ | omb_ |
obs - mean(H(x)); set in computeHofX method More... | |
DeparturesEnsemble_ | Yb_ |
std::unique_ptr< ObsErrors_ > | R_ |
observation errors, set in computeHofX method More... | |
std::unique_ptr< Departures_ > | invVarR_ |
Private Types | |
typedef Departures< OBS > | Departures_ |
typedef DeparturesEnsemble< OBS > | DeparturesEnsemble_ |
typedef Geometry< MODEL > | Geometry_ |
typedef GeometryIterator< MODEL > | GeometryIterator_ |
typedef IncrementEnsemble4D< MODEL > | IncrementEnsemble4D_ |
typedef ObsErrors< OBS > | ObsErrors_ |
typedef ObsLocalizations< MODEL, OBS > | ObsLocalizations_ |
typedef ObsSpaces< OBS > | ObsSpaces_ |
Local Ensemble Tranform Kalman Filter solver.
An implementation of the LETKF from Hunt et al. 2007 this version is implemented using Eigen algebra and temporary Eigen matrices for Xa and Xb this verion implements RTPP and RTPS.
Hunt, B. R., Kostelich, E. J., & Szunyogh, I. (2007). Efficient data assimilation for spatiotemporal chaos: A local ensemble transform Kalman filter. Physica D: Nonlinear Phenomena, 230(1-2), 112-126.
Definition at line 47 of file LETKFSolver.h.
|
private |
Definition at line 48 of file LETKFSolver.h.
|
private |
Definition at line 49 of file LETKFSolver.h.
|
private |
Definition at line 50 of file LETKFSolver.h.
|
private |
Definition at line 51 of file LETKFSolver.h.
|
private |
Definition at line 52 of file LETKFSolver.h.
|
private |
Definition at line 53 of file LETKFSolver.h.
|
private |
Definition at line 54 of file LETKFSolver.h.
|
private |
Definition at line 55 of file LETKFSolver.h.
oops::LETKFSolver< MODEL, OBS >::LETKFSolver | ( | ObsSpaces_ & | obspaces, |
const Geometry_ & | geometry, | ||
const eckit::Configuration & | config, | ||
size_t | nens | ||
) |
|
protectedvirtual |
Applies weights and adds posterior inflation.
Definition at line 202 of file LETKFSolver.h.
|
inlinestatic |
Definition at line 58 of file LETKFSolver.h.
|
protectedvirtual |
Computes weights for ensemble update with local observations
[in] | omb | Observation departures (nlocalobs) |
[in] | Yb | Ensemble perturbations (nens, nlocalobs) |
[in] | invvarR | Inverse of observation error variances (nlocalobs) |
Reimplemented in oops::LETKFSolverGSI< MODEL, OBS >.
Definition at line 166 of file LETKFSolver.h.
|
overridevirtual |
KF update + posterior inflation at a grid point location (GeometryIterator_)
Implements oops::LocalEnsembleSolver< MODEL, OBS >.
Definition at line 133 of file LETKFSolver.h.
|
protected |
Definition at line 84 of file LETKFSolver.h.
|
protected |
Definition at line 85 of file LETKFSolver.h.
|
protected |
Definition at line 87 of file LETKFSolver.h.
|
protected |
Definition at line 78 of file LETKFSolver.h.
|
protected |
Definition at line 80 of file LETKFSolver.h.
|
protected |
Definition at line 81 of file LETKFSolver.h.