#include <KalmanSmoother.h>
Implements the Rauch–Tung–Striebel smoother.
Definition at line 14 of file KalmanSmoother.h.
|
static std::pair< Eigen::VectorXd, Eigen::MatrixXd > | Analyze (std::pair< Eigen::VectorXd, Eigen::MatrixXd > const &currDist_t, std::pair< Eigen::VectorXd, Eigen::MatrixXd > const &nextDist_t, std::pair< Eigen::VectorXd, Eigen::MatrixXd > const &nextDist_n, std::shared_ptr< muq::Modeling::LinearOperator > F) |
|
template<typename MatrixType > |
static std::pair< Eigen::VectorXd, Eigen::MatrixXd > | Analyze (std::pair< Eigen::VectorXd, Eigen::MatrixXd > const &currDist_t, std::pair< Eigen::VectorXd, Eigen::MatrixXd > const &nextDist_t, std::pair< Eigen::VectorXd, Eigen::MatrixXd > const &nextDist_n, MatrixType const &F) |
|
◆ Analyze() [1/2]
template<typename MatrixType >
static std::pair<Eigen::VectorXd, Eigen::MatrixXd> muq::Inference::KalmanSmoother::Analyze |
( |
std::pair< Eigen::VectorXd, Eigen::MatrixXd > const & |
currDist_t, |
|
|
std::pair< Eigen::VectorXd, Eigen::MatrixXd > const & |
nextDist_t, |
|
|
std::pair< Eigen::VectorXd, Eigen::MatrixXd > const & |
nextDist_n, |
|
|
MatrixType const & |
F |
|
) |
| |
|
inlinestatic |
- Template Parameters
-
MatrixType | A type that can be converted to a MUQ LinearOperator. Examples include Eigen::MatrixXd and Eigen::SparseMatrix |
Definition at line 34 of file KalmanSmoother.h.
References Analyze().
◆ Analyze() [2/2]
std::pair< Eigen::VectorXd, Eigen::MatrixXd > KalmanSmoother::Analyze |
( |
std::pair< Eigen::VectorXd, Eigen::MatrixXd > const & |
currDist_t, |
|
|
std::pair< Eigen::VectorXd, Eigen::MatrixXd > const & |
nextDist_t, |
|
|
std::pair< Eigen::VectorXd, Eigen::MatrixXd > const & |
nextDist_n, |
|
|
std::shared_ptr< muq::Modeling::LinearOperator > |
F |
|
) |
| |
|
static |
- Parameters
-
[in] | currDist_t | The distribution at time t after the forward Kalman filtering step (i.e., using all data up to and including time t). |
[in] | nextDist_t | The distribution at time t+1 from the prediction phase of the Kalman filtering step (i.e., using all data up to time t). |
[in] | nextDist_n | The distribution at time t+1 after the RTS smoothing step (i.e., using all data). |
[in] | F | The linear operator acting on the state at time t, to produce the state at time t+1 |
- Returns
- A distribution (mean and covariance) at time t after accounting for all data.
Definition at line 9 of file KalmanSmoother.cpp.
Referenced by Analyze(), and muq::Inference::PythonBindings::KalmanWrapper().
◆ ComputeC()
Eigen::MatrixXd KalmanSmoother::ComputeC |
( |
Eigen::MatrixXd const & |
currDist_t_cov, |
|
|
Eigen::MatrixXd const & |
nextDist_t_cov, |
|
|
std::shared_ptr< muq::Modeling::LinearOperator > const & |
F |
|
) |
| |
|
staticprivate |
The documentation for this class was generated from the following files: