00001 #ifndef INDII_ML_FILTER_UNSCENTEDKALMANSMOOTHER_HPP
00002 #define INDII_ML_FILTER_UNSCENTEDKALMANSMOOTHER_HPP
00003
00004 #include "UnscentedKalmanFilter.hpp"
00005 #include "TwoFilterSmoother.hpp"
00006
00007 namespace indii {
00008 namespace ml {
00009 namespace filter {
00010
00011 template<class T> class UnscentedKalmanSmootherModel;
00012 template<class T> class UnscentedKalmanSmootherBackwardTransitionAdaptor;
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033 template <class T>
00034 class UnscentedKalmanSmoother : public TwoFilterSmoother<T> {
00035
00036 friend class UnscentedKalmanSmootherBackwardTransitionAdaptor<T>;
00037
00038 public:
00039
00040
00041
00042
00043
00044
00045
00046
00047 UnscentedKalmanSmoother(UnscentedKalmanSmootherModel<T>* model,
00048 const T tT,
00049 const indii::ml::aux::GaussianPdf& p_xT);
00050
00051
00052
00053
00054 virtual ~UnscentedKalmanSmoother();
00055
00056 virtual void smooth(const T tn, const indii::ml::aux::vector& ytn,
00057 const indii::ml::aux::GaussianPdf& p_xtn_ytn);
00058
00059 virtual indii::ml::aux::GaussianPdf backwardMeasure();
00060
00061 virtual indii::ml::aux::GaussianPdf smoothedMeasure();
00062
00063
00064
00065
00066
00067 UnscentedTransformation<T>* getBackwardTransitionTransformation();
00068
00069
00070
00071
00072
00073 UnscentedTransformation<T>* getMeasurementTransformation();
00074
00075 protected:
00076
00077
00078
00079 const unsigned int N;
00080
00081
00082
00083
00084 const unsigned int W;
00085
00086
00087
00088
00089 const unsigned int V;
00090
00091
00092
00093
00094 aux::GaussianPdf p_Xtn_ytn_b;
00095
00096
00097
00098
00099 UnscentedKalmanSmootherModel<T>* model;
00100
00101
00102
00103
00104 UnscentedTransformationModel<T>* backwardTransitionModel;
00105
00106
00107
00108
00109 UnscentedTransformation<T>* backwardTransitionTransform;
00110
00111
00112
00113
00114 UnscentedTransformationModel<T>* measurementModel;
00115
00116
00117
00118
00119 UnscentedTransformation<T>* measurementTransform;
00120
00121 };
00122
00123 }
00124 }
00125 }
00126
00127 #include "UnscentedKalmanSmootherBackwardTransitionAdaptor.hpp"
00128
00129 #include "boost/numeric/bindings/traits/ublas_matrix.hpp"
00130 #include "boost/numeric/bindings/traits/ublas_vector.hpp"
00131 #include "boost/numeric/bindings/lapack/lapack.hpp"
00132
00133 template <class T>
00134 indii::ml::filter::UnscentedKalmanSmoother<T>::UnscentedKalmanSmoother(
00135 UnscentedKalmanSmootherModel<T>* model, const T tT,
00136 const aux::GaussianPdf& p_xT) :
00137 TwoFilterSmoother<T>(tT, p_xT),
00138 N(model->getStateSize()),
00139 W(model->getSystemNoise().getDimensions()),
00140 V(model->getMeasurementNoise().getDimensions()),
00141 p_Xtn_ytn_b(N+W+V),
00142 model(model) {
00143 namespace aux = indii::ml::aux;
00144 namespace ublas = boost::numeric::ublas;
00145
00146
00147 assert (p_xT.getDimensions() == model->getStateSize());
00148
00149
00150 backwardTransitionModel =
00151 new UnscentedKalmanSmootherBackwardTransitionAdaptor<T>(model);
00152 backwardTransitionTransform = new UnscentedTransformation<T>(
00153 *backwardTransitionModel);
00154 measurementModel = new UnscentedKalmanFilterMeasurementAdaptor<T>(model);
00155 measurementTransform = new UnscentedTransformation<T>(*measurementModel);
00156
00157
00158 ublas::range xRange(0,N);
00159 ublas::range wRange(N,N+W);
00160 ublas::range vRange(N+W,N+W+V);
00161
00162
00163 aux::vector mu(N+W+V);
00164 aux::symmetric_matrix sigma(N+W+V);
00165 aux::GaussianPdf& w = model->getSystemNoise();
00166 aux::GaussianPdf& v = model->getMeasurementNoise();
00167
00168 noalias(project(mu, xRange)) = p_xT.getExpectation();
00169 noalias(project(mu, wRange)) = w.getExpectation();
00170 noalias(project(mu, vRange)) = v.getExpectation();
00171
00172 noalias(project(sigma, xRange, xRange)) = p_xT.getCovariance();
00173 noalias(project(sigma, wRange, xRange)) = aux::zero_matrix(W, N);
00174 noalias(project(sigma, wRange, wRange)) = w.getCovariance();
00175 noalias(project(sigma, vRange, xRange)) = aux::zero_matrix(V, N);
00176 noalias(project(sigma, vRange, wRange)) = aux::zero_matrix(V, W);
00177 noalias(project(sigma, vRange, vRange)) = v.getCovariance();
00178
00179 p_Xtn_ytn_b.setExpectation(mu);
00180 p_Xtn_ytn_b.setCovariance(sigma);
00181 }
00182
00183 template <class T>
00184 indii::ml::filter::UnscentedKalmanSmoother<T>::~UnscentedKalmanSmoother() {
00185 delete backwardTransitionTransform;
00186 delete backwardTransitionModel;
00187 delete measurementTransform;
00188 delete measurementModel;
00189 }
00190
00191 template <class T>
00192 void indii::ml::filter::UnscentedKalmanSmoother<T>::smooth(const T tn,
00193 const indii::ml::aux::vector& ytn,
00194 const indii::ml::aux::GaussianPdf& p_xtn_ytn) {
00195 namespace aux = indii::ml::aux;
00196 namespace lapack = boost::numeric::bindings::lapack;
00197 namespace ublas = boost::numeric::ublas;
00198
00199
00200 assert (tn < this->tn);
00201
00202 T delta = this->tn - tn;
00203 this->tn = tn;
00204
00205
00206 ublas::range xRange(0,N);
00207 ublas::range wRange(N,N+W);
00208 ublas::range vRange(N+W,N+W+V);
00209
00210
00211
00212 aux::GaussianPdf p_xtn_ytnp1_b(backwardTransitionTransform->transform(
00213 p_Xtn_ytn_b, delta));
00214
00215
00216 aux::vector mu(p_Xtn_ytn_b.getExpectation());
00217 aux::symmetric_matrix sigma(p_Xtn_ytn_b.getCovariance());
00218
00219 noalias(project(mu, xRange)) = p_xtn_ytnp1_b.getExpectation();
00220 noalias(project(sigma, xRange, xRange)) = p_xtn_ytnp1_b.getCovariance();
00221
00222 aux::GaussianPdf p_Xtn_ytnp1(mu, sigma);
00223
00224
00225
00226 aux::matrix P_xy(N+W+V,V);
00227 aux::GaussianPdf p_ytn_Xtn(
00228 measurementTransform->transform(p_Xtn_ytnp1, delta, &P_xy));
00229
00230
00231 const aux::vector &y_ytn_Xtn = p_ytn_Xtn.getExpectation();
00232 const aux::symmetric_matrix &P_ytn_Xtn = p_ytn_Xtn.getCovariance();
00233 const aux::vector &x_Xtn_ytnp1 = p_Xtn_ytnp1.getExpectation();
00234 const aux::symmetric_matrix &P_Xtn_ytnp1 = p_Xtn_ytnp1.getCovariance();
00235 aux::matrix PI_ytn_Xtn(V,V);
00236 aux::matrix X(P_ytn_Xtn);
00237 aux::inv(X, PI_ytn_Xtn);
00238 aux::matrix K(prod(P_xy, PI_ytn_Xtn));
00239
00240
00241
00242 noalias(project(K,wRange,ublas::range(0,V))) = aux::zero_matrix(W,V);
00243 noalias(project(K,vRange,ublas::range(0,V))) = aux::zero_matrix(V,V);
00244
00245 noalias(mu) = x_Xtn_ytnp1 + prod(K, ytn - y_ytn_Xtn);
00246 X.resize(V, N+W+V, false);
00247 noalias(X) = prod(P_ytn_Xtn, trans(K));
00248 noalias(sigma) = P_Xtn_ytnp1 - prod(K, X);
00249
00250 p_Xtn_ytn_b.setExpectation(mu);
00251 p_Xtn_ytn_b.setCovariance(sigma);
00252 this->p_xtn_ytn_b.setExpectation(project(mu, xRange));
00253 this->p_xtn_ytn_b.setCovariance(project(sigma, xRange,
00254 xRange));
00255
00256
00257 const aux::symmetric_matrix& P_xtn_ytn= p_xtn_ytn.getCovariance();
00258 const aux::symmetric_matrix& P_xtn_ytnp1_b = p_xtn_ytnp1_b.getCovariance();
00259 aux::matrix PI_xtn_ytn(N,N);
00260 aux::matrix PI_xtn_ytnp1_b(N,N);
00261 X.resize(N,N,false);
00262
00263
00264 noalias(X) = P_xtn_ytn;
00265 aux::inv(X, PI_xtn_ytn);
00266
00267 noalias(X) = P_xtn_ytnp1_b;
00268 aux::inv(X, PI_xtn_ytnp1_b);
00269
00270
00271 aux::symmetric_matrix P_xtn_ytT(N);
00272 aux::matrix Y(N,N);
00273 noalias(Y) = PI_xtn_ytn + PI_xtn_ytnp1_b;
00274
00275 aux::inv(Y,X);
00276 noalias(P_xtn_ytT) = ublas::symmetric_adaptor<aux::matrix, ublas::lower>(X);
00277
00278
00279 const aux::vector& x_xtn_ytn = p_xtn_ytn.getExpectation();
00280 const aux::vector& x_xtn_ytnp1_b = p_xtn_ytnp1_b.getExpectation();
00281 aux::vector x_xtn_ytT(N);
00282
00283 noalias(x_xtn_ytT) = prod(PI_xtn_ytn, x_xtn_ytn);
00284 noalias(x_xtn_ytT) += prod(PI_xtn_ytnp1_b, x_xtn_ytnp1_b);
00285 x_xtn_ytT = prod(P_xtn_ytT, x_xtn_ytT);
00286
00287
00288 this->p_xtn_ytT.setExpectation(x_xtn_ytT);
00289 this->p_xtn_ytT.setCovariance(P_xtn_ytT);
00290
00291
00292 assert (this->tn == tn);
00293 }
00294
00295 template <class T>
00296 indii::ml::aux::GaussianPdf
00297 indii::ml::filter::UnscentedKalmanSmoother<T>::backwardMeasure() {
00298 return measurementTransform->transform(p_Xtn_ytn_b);
00299 }
00300
00301 template <class T>
00302 indii::ml::aux::GaussianPdf
00303 indii::ml::filter::UnscentedKalmanSmoother<T>::smoothedMeasure() {
00304 namespace aux = indii::ml::aux;
00305 namespace ublas = boost::numeric::ublas;
00306
00307
00308 ublas::range xRange(0,N);
00309 ublas::range wRange(N,N+W);
00310 ublas::range vRange(N+W,N+W+V);
00311
00312
00313 aux::vector mu(p_Xtn_ytn_b.getExpectation());
00314 aux::symmetric_matrix sigma(p_Xtn_ytn_b.getCovariance());
00315
00316 noalias(project(mu, xRange)) = this->p_xtn_ytT.getExpectation();
00317 noalias(project(sigma, xRange, xRange)) =
00318 this->p_xtn_ytT.getCovariance();
00319
00320 aux::GaussianPdf p_Xtn_ytT(mu, sigma);
00321
00322 return measurementTransform->transform(p_Xtn_ytT);
00323 }
00324
00325 template <class T>
00326 indii::ml::filter::UnscentedTransformation<T>*
00327 indii::ml::filter::UnscentedKalmanSmoother<T>::getBackwardTransitionTransformation() {
00328 return backwardTransitionTransform;
00329 }
00330
00331 template <class T>
00332 indii::ml::filter::UnscentedTransformation<T>*
00333 indii::ml::filter::UnscentedKalmanSmoother<T>::getMeasurementTransformation() {
00334 return measurementTransform;
00335 }
00336
00337 #endif
00338