Bayes Filters Library
SIS.cpp
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2016-2019 Istituto Italiano di Tecnologia (IIT)
3  *
4  * This software may be modified and distributed under the terms of the
5  * BSD 3-Clause license. See the accompanying LICENSE file for details.
6  */
7 
8 #include <BayesFilters/SIS.h>
9 #include <BayesFilters/utils.h>
10 
11 #include <fstream>
12 #include <iostream>
13 #include <utility>
14 
15 #include <Eigen/Dense>
16 
17 using namespace bfl;
18 using namespace Eigen;
19 
20 
22 (
23  unsigned int num_particle,
24  std::size_t state_size_linear,
25  std::size_t state_size_circular,
26  std::unique_ptr<ParticleSetInitialization> initialization,
27  std::unique_ptr<PFPrediction> prediction,
28  std::unique_ptr<PFCorrection> correction,
29  std::unique_ptr<Resampling> resampling
30 ) noexcept :
31  ParticleFilter(std::move(initialization), std::move(prediction), std::move(correction), std::move(resampling)),
32  num_particle_(num_particle),
33  state_size_(state_size_linear + state_size_circular),
34  pred_particle_(num_particle_, state_size_linear, state_size_circular),
35  cor_particle_(num_particle_, state_size_linear, state_size_circular)
36 { }
37 
38 
40 (
41  unsigned int num_particle,
42  std::size_t state_size_linear,
43  std::unique_ptr<ParticleSetInitialization> initialization,
44  std::unique_ptr<PFPrediction> prediction,
45  std::unique_ptr<PFCorrection> correction,
46  std::unique_ptr<Resampling> resampling
47 ) noexcept :
48  SIS(num_particle, state_size_linear, 0, std::move(initialization), std::move(prediction), std::move(correction), std::move(resampling))
49 { }
50 
51 
53 {
54  return initialization().initialize(pred_particle_);
55 }
56 
57 
59 {
60  if (step_number() != 0)
61  prediction().predict(cor_particle_, pred_particle_);
62 
63  if (correction().freeze_measurements())
64  {
65  correction().correct(pred_particle_, cor_particle_);
66 
67  /* Normalize weights using LogSumExp. */
68  cor_particle_.weight().array() -= utils::log_sum_exp(cor_particle_.weight());
69  }
70  else
71  cor_particle_ = pred_particle_;
72 
73  log();
74 
75  if (resampling().neff(cor_particle_.weight()) < static_cast<double>(num_particle_)/3.0)
76  {
77  ParticleSet res_particle(num_particle_, state_size_);
78  VectorXi res_parent(num_particle_, 1);
79 
80  resampling().resample(cor_particle_, res_particle, res_parent);
81 
82  cor_particle_ = res_particle;
83  }
84 }
85 
86 
88 {
89  return true;
90 }
91 
92 
93 std::vector<std::string> SIS::log_file_names(const std::string& folder_path, const std::string& file_name_prefix)
94 {
95  return { folder_path + "/" + file_name_prefix + "_pred_particles",
96  folder_path + "/" + file_name_prefix + "_pred_weights",
97  folder_path + "/" + file_name_prefix + "_cor_particles",
98  folder_path + "/" + file_name_prefix + "_cor_weights" };
99 }
100 
101 
102 void SIS::log()
103 {
104  logger(pred_particle_.state().transpose(), pred_particle_.weight().transpose(),
105  cor_particle_.state().transpose(), cor_particle_.weight().transpose());
106 }
bfl::utils::log_sum_exp
double log_sum_exp(const Eigen::MatrixBase< Derived > &data)
Return the element-wise logarithm of the sum of exponentials of the input data.
Definition: utils.h:72
bfl::SIS::log
void log() override
Definition: SIS.cpp:102
bfl::SIS
Definition: SIS.h:25
bfl::SIS::run_condition
bool run_condition() override
Definition: SIS.cpp:87
bfl
Port of boost::any for C++11 compilers.
Definition: AdditiveMeasurementModel.h:13
bfl::SIS::initialization_step
bool initialization_step() override
Definition: SIS.cpp:52
bfl::SIS::log_file_names
std::vector< std::string > log_file_names(const std::string &folder_path, const std::string &file_name_prefix) override
Definition: SIS.cpp:93
SIS.h
utils.h
bfl::SIS::SIS
SIS(unsigned int num_particle, std::size_t state_size_linear, std::size_t state_size_circular, std::unique_ptr< ParticleSetInitialization > initialization, std::unique_ptr< PFPrediction > prediction, std::unique_ptr< PFCorrection > correction, std::unique_ptr< Resampling > resampling) noexcept
Definition: SIS.cpp:22
bfl::SIS::filtering_step
void filtering_step() override
Definition: SIS.cpp:58
bfl::ParticleFilter
Definition: ParticleFilter.h:24
bfl::ParticleSet
Definition: ParticleSet.h:20