LArSoft  v10_06_00
Liquid Argon Software toolkit - https://larsoft.org/
lar_content::KalmanFilter3D Class Reference

KalmanFilter3D class. More...

#include "KalmanFilter.h"

Public Member Functions

 KalmanFilter3D (double dt, double processVariance, double measurementVariance, const Eigen::VectorXd &x)
 Default constructor. More...
 
 ~KalmanFilter3D ()
 Destructor. More...
 
void Predict ()
 Uses the state transition matrix to predict the next state. More...
 
void Update (const Eigen::VectorXd &z)
 Uses the measurement matrix to update the state estimate. More...
 
const Eigen::VectorXd & GetState () const
 Get the current state vector. More...
 
const Eigen::VectorXd GetPosition () const
 Get the current position vector. More...
 
const Eigen::VectorXd GetDirection () const
 Get the current direction vector. More...
 
const Eigen::VectorXd & GetTemporaryState () const
 Get the temporary state vector. More...
 

Private Attributes

double m_dt
 Time step. More...
 
Eigen::VectorXd m_x
 State vector (tracks position and 'velocity') More...
 
Eigen::MatrixXd m_P
 Covariance matrix. More...
 
Eigen::MatrixXd m_F
 State transition matrix. More...
 
Eigen::MatrixXd m_H
 Measurement matrix. More...
 
Eigen::MatrixXd m_R
 Measurement covariance matrix. More...
 
Eigen::MatrixXd m_Q
 Process covariance matrix. More...
 
Eigen::VectorXd m_xTemp
 Temporary state vector. More...
 
Eigen::MatrixXd m_PTemp
 Temporary covariance matrix. More...
 

Detailed Description

KalmanFilter3D class.

Definition at line 19 of file KalmanFilter.h.

Constructor & Destructor Documentation

lar_content::KalmanFilter3D::KalmanFilter3D ( double  dt,
double  processVariance,
double  measurementVariance,
const Eigen::VectorXd &  x 
)

Default constructor.

Parameters
dtTime step
processVarianceProcess variance
measurementVarianceMeasurement variance
xInitial position

Definition at line 16 of file KalmanFilter.cc.

References m_F, m_H, m_P, m_Q, m_R, m_x, and x.

16  :
17  m_dt{dt}
18 {
19  m_x = VectorXd::Zero(6);
20  m_x.head(3) = x;
21  m_P = MatrixXd::Identity(6, 6);
22  m_F = MatrixXd::Identity(6, 6);
23  for (int i = 0; i < 3; ++i)
24  {
25  m_F(i, i + 3) = dt;
26  }
27  m_H = MatrixXd::Zero(3, 6);
28  m_H.block<3, 3>(0, 0) = MatrixXd::Identity(3, 3);
29  m_Q = MatrixXd::Zero(6, 6);
30  for (int i = 0; i < 3; ++i)
31  {
32  m_Q(i, i) = (0.25 * dt * dt * dt * dt) * processVariance;
33  m_Q(i, i + 3) = (0.5 * dt * dt * dt) * processVariance;
34  m_Q(i + 3, i) = (0.5 * dt * dt * dt) * processVariance;
35  m_Q(i + 3, i + 3) = (dt * dt) * processVariance;
36  }
37  m_R = MatrixXd::Identity(3, 3) * measurementVariance;
38 }
Float_t x
Definition: compare.C:6
Eigen::MatrixXd m_R
Measurement covariance matrix.
Definition: KalmanFilter.h:83
Eigen::MatrixXd m_P
Covariance matrix.
Definition: KalmanFilter.h:80
Eigen::MatrixXd m_H
Measurement matrix.
Definition: KalmanFilter.h:82
Eigen::VectorXd m_x
State vector (tracks position and &#39;velocity&#39;)
Definition: KalmanFilter.h:79
Eigen::MatrixXd m_F
State transition matrix.
Definition: KalmanFilter.h:81
double m_dt
Time step.
Definition: KalmanFilter.h:78
Eigen::MatrixXd m_Q
Process covariance matrix.
Definition: KalmanFilter.h:84
lar_content::KalmanFilter3D::~KalmanFilter3D ( )

Destructor.

Definition at line 42 of file KalmanFilter.cc.

43 {
44 }

Member Function Documentation

const VectorXd lar_content::KalmanFilter3D::GetDirection ( ) const

Get the current direction vector.

Returns
The current direction vector

Definition at line 83 of file KalmanFilter.cc.

References m_x.

84 {
85  return m_x.tail(3).normalized();
86 }
Eigen::VectorXd m_x
State vector (tracks position and &#39;velocity&#39;)
Definition: KalmanFilter.h:79
const VectorXd lar_content::KalmanFilter3D::GetPosition ( ) const

Get the current position vector.

Returns
The current postion vector

Definition at line 76 of file KalmanFilter.cc.

References m_x.

77 {
78  return m_x.head(3);
79 }
Eigen::VectorXd m_x
State vector (tracks position and &#39;velocity&#39;)
Definition: KalmanFilter.h:79
const VectorXd & lar_content::KalmanFilter3D::GetState ( ) const

Get the current state vector.

Returns
The current state vector

Definition at line 69 of file KalmanFilter.cc.

References m_x.

70 {
71  return m_x;
72 }
Eigen::VectorXd m_x
State vector (tracks position and &#39;velocity&#39;)
Definition: KalmanFilter.h:79
const VectorXd & lar_content::KalmanFilter3D::GetTemporaryState ( ) const

Get the temporary state vector.

Returns
The temporary state vector

Definition at line 90 of file KalmanFilter.cc.

References m_xTemp.

91 {
92  return m_xTemp;
93 }
Eigen::VectorXd m_xTemp
Temporary state vector.
Definition: KalmanFilter.h:85
void lar_content::KalmanFilter3D::Predict ( )

Uses the state transition matrix to predict the next state.

Definition at line 48 of file KalmanFilter.cc.

References m_F, m_P, m_PTemp, m_Q, m_x, and m_xTemp.

49 {
50  m_xTemp = m_F * m_x;
51  m_PTemp = m_F * m_P * m_F.transpose() + m_Q;
52 }
Eigen::MatrixXd m_P
Covariance matrix.
Definition: KalmanFilter.h:80
Eigen::VectorXd m_xTemp
Temporary state vector.
Definition: KalmanFilter.h:85
Eigen::VectorXd m_x
State vector (tracks position and &#39;velocity&#39;)
Definition: KalmanFilter.h:79
Eigen::MatrixXd m_PTemp
Temporary covariance matrix.
Definition: KalmanFilter.h:86
Eigen::MatrixXd m_F
State transition matrix.
Definition: KalmanFilter.h:81
Eigen::MatrixXd m_Q
Process covariance matrix.
Definition: KalmanFilter.h:84
void lar_content::KalmanFilter3D::Update ( const Eigen::VectorXd &  z)

Uses the measurement matrix to update the state estimate.

Parameters
zThe position to use for updating the state estimate

Definition at line 56 of file KalmanFilter.cc.

References m_H, m_P, m_PTemp, m_R, m_x, m_xTemp, and y.

57 {
58  m_x = m_xTemp;
59  m_P = m_PTemp;
60  VectorXd y = z - m_H * m_x;
61  MatrixXd S = m_H * m_P * m_H.transpose() + m_R;
62  MatrixXd K = m_P * m_H.transpose() * S.inverse();
63  m_x += K * y;
64  m_P = (MatrixXd::Identity(6, 6) - K * m_H) * m_P;
65 }
Eigen::MatrixXd m_R
Measurement covariance matrix.
Definition: KalmanFilter.h:83
Eigen::MatrixXd m_P
Covariance matrix.
Definition: KalmanFilter.h:80
Float_t y
Definition: compare.C:6
Double_t z
Definition: plot.C:276
Eigen::MatrixXd m_H
Measurement matrix.
Definition: KalmanFilter.h:82
Eigen::VectorXd m_xTemp
Temporary state vector.
Definition: KalmanFilter.h:85
Eigen::VectorXd m_x
State vector (tracks position and &#39;velocity&#39;)
Definition: KalmanFilter.h:79
Eigen::MatrixXd m_PTemp
Temporary covariance matrix.
Definition: KalmanFilter.h:86

Member Data Documentation

double lar_content::KalmanFilter3D::m_dt
private

Time step.

Definition at line 78 of file KalmanFilter.h.

Eigen::MatrixXd lar_content::KalmanFilter3D::m_F
private

State transition matrix.

Definition at line 81 of file KalmanFilter.h.

Referenced by KalmanFilter3D(), and Predict().

Eigen::MatrixXd lar_content::KalmanFilter3D::m_H
private

Measurement matrix.

Definition at line 82 of file KalmanFilter.h.

Referenced by KalmanFilter3D(), and Update().

Eigen::MatrixXd lar_content::KalmanFilter3D::m_P
private

Covariance matrix.

Definition at line 80 of file KalmanFilter.h.

Referenced by KalmanFilter3D(), Predict(), and Update().

Eigen::MatrixXd lar_content::KalmanFilter3D::m_PTemp
private

Temporary covariance matrix.

Definition at line 86 of file KalmanFilter.h.

Referenced by Predict(), and Update().

Eigen::MatrixXd lar_content::KalmanFilter3D::m_Q
private

Process covariance matrix.

Definition at line 84 of file KalmanFilter.h.

Referenced by KalmanFilter3D(), and Predict().

Eigen::MatrixXd lar_content::KalmanFilter3D::m_R
private

Measurement covariance matrix.

Definition at line 83 of file KalmanFilter.h.

Referenced by KalmanFilter3D(), and Update().

Eigen::VectorXd lar_content::KalmanFilter3D::m_x
private

State vector (tracks position and 'velocity')

Definition at line 79 of file KalmanFilter.h.

Referenced by GetDirection(), GetPosition(), GetState(), KalmanFilter3D(), Predict(), and Update().

Eigen::VectorXd lar_content::KalmanFilter3D::m_xTemp
private

Temporary state vector.

Definition at line 85 of file KalmanFilter.h.

Referenced by GetTemporaryState(), Predict(), and Update().


The documentation for this class was generated from the following files: