Belle II Software development
KalmanCalculator Class Reference

does the calculation of the gain matrix, updates the cov and fitpars More...

#include <KalmanCalculator.h>

Public Member Functions

 KalmanCalculator (int sizeRes, int sizeState)
 constructor
 
void updateState (FitParams &fitparams)
 update statevector
 
void updateState (FitParams &fitparams, FitParams &oldState)
 update statevector
 
void updateCovariance (FitParams &fitparams)
 update the statevectors covariance
 
double getChiSquare () const
 get chi2 of this iteration
 
ErrCode calculateGainMatrix (const Eigen::Matrix< double, -1, 1, 0, 7, 1 > &residuals, const Eigen::Matrix< double, -1, -1, 0, 7, MAX_MATRIX_SIZE > &G, const FitParams &fitparams, const Eigen::Matrix< double, -1, -1, 0, 7, 7 > *V=0, double weight=1)
 init the kalman machienery
 
double getConstraintDim () const
 get dimension of the constraint
 

Private Attributes

int m_constrDim
 dimension of the constraint
 
double m_chisq
 chi2
 
Eigen::Matrix< double, -1, 1, 0, 7, 1 > m_res
 we know the max sizes of the matrices we assume the tree is smaller than MAX_MATRIX_SIZE parameters and the largest constraint is the track constraint with 7
 
Eigen::Matrix< double, -1, -1, 0, 7, MAX_MATRIX_SIZE > m_G
 G former H, transforms covraince of {residuals}<->{x,p,E}.
 
Eigen::Matrix< double, -1, -1, 0, 7, 7 > m_R
 R residual covariance.
 
Eigen::Matrix< double, -1, -1, 0, 7, 7 > m_Rinverse
 R inverse.
 
Eigen::Matrix< double, -1, -1, 0, MAX_MATRIX_SIZE, 7 > m_K
 K kalman gain matrix.
 
Eigen::Matrix< double, -1, -1, 0, MAX_MATRIX_SIZE, 7 > m_CGt
 C times G^t.
 

Detailed Description

does the calculation of the gain matrix, updates the cov and fitpars

Definition at line 22 of file KalmanCalculator.h.

Constructor & Destructor Documentation

◆ KalmanCalculator()

KalmanCalculator ( int sizeRes,
int sizeState )

constructor

Definition at line 15 of file KalmanCalculator.cc.

18 :
19 m_constrDim(sizeRes),
20 m_chisq(1e10),
21 m_res(sizeRes),
22 m_G(sizeRes, sizeState),
23 m_R(sizeRes, sizeRes),
24 m_Rinverse(sizeRes, sizeRes),
25 m_K(sizeState, sizeRes),
26 m_CGt(sizeState, sizeRes)
27 {
28 m_R = Eigen::Matrix < double, -1, -1, 0, 7, 7 >::Zero(m_constrDim, m_constrDim);
29 }

Member Function Documentation

◆ calculateGainMatrix()

ErrCode calculateGainMatrix ( const Eigen::Matrix< double, -1, 1, 0, 7, 1 > & residuals,
const Eigen::Matrix< double, -1, -1, 0, 7, MAX_MATRIX_SIZE > & G,
const FitParams & fitparams,
const Eigen::Matrix< double, -1, -1, 0, 7, 7 > * V = 0,
double weight = 1 )

init the kalman machienery

Definition at line 32 of file KalmanCalculator.cc.

38 {
39 m_res = residuals;
40 m_G = G;
41
42 Eigen::Matrix < double, -1, -1, 0, MAX_MATRIX_SIZE, MAX_MATRIX_SIZE > C = fitparams.getCovariance().triangularView<Eigen::Lower>();
43
44 m_CGt = C.selfadjointView<Eigen::Lower>() * G.transpose();
45 Eigen::Matrix < double, -1, -1, 0, 7, 7 > Rtemp = G * m_CGt;
46 if (V && (weight) && ((*V).diagonal().array() != 0).all()) {
47
48 const Eigen::Matrix < double, -1, -1, 0, 7, 7 > weightedV =
49 weight * (*V).selfadjointView<Eigen::Lower>();
50
51 m_R = Rtemp + weightedV;
52
53 } else {
54 m_R = Rtemp.triangularView<Eigen::Lower>();
55 }
56
57 Eigen::Matrix < double, -1, -1, 0, 7, 7 > RInvtemp;
58 RInvtemp = m_R.selfadjointView<Eigen::Lower>();
59 m_Rinverse = RInvtemp.inverse();
60 if (!m_Rinverse.allFinite()) { return ErrCode(ErrCode::Status::inversionerror); }
61
62 m_K = m_CGt * m_Rinverse.selfadjointView<Eigen::Lower>();
63 return ErrCode(ErrCode::Status::success);
64 }

◆ getChiSquare()

double getChiSquare ( ) const
inline

get chi2 of this iteration

Definition at line 41 of file KalmanCalculator.h.

41{ return m_chisq;}

◆ getConstraintDim()

double getConstraintDim ( ) const
inline

get dimension of the constraint

Definition at line 52 of file KalmanCalculator.h.

52{ return m_constrDim; }

◆ updateCovariance()

TREEFITTER_NO_STACK_WARNING void updateCovariance ( FitParams & fitparams)

update the statevectors covariance

Definition at line 88 of file KalmanCalculator.cc.

89 {
90 Eigen::Matrix < double, -1, -1, 0, MAX_MATRIX_SIZE, MAX_MATRIX_SIZE > fitCov =
91 fitparams.getCovariance().triangularView<Eigen::Lower>();
92
93 Eigen::Matrix < double, -1, -1, 0, MAX_MATRIX_SIZE, MAX_MATRIX_SIZE > GRinvGt =
94 m_G.transpose() * m_Rinverse.selfadjointView<Eigen::Lower>() * m_G;
95
96 //fitcov is sym so no transpose needed (not that it would have runtime cost)
97 Eigen::Matrix < double, -1, -1, 0, MAX_MATRIX_SIZE, MAX_MATRIX_SIZE > deltaCov =
98 fitCov.selfadjointView<Eigen::Lower>() * GRinvGt * fitCov.selfadjointView<Eigen::Lower>();
99
100 double eps = Eigen::NumTraits<double>::epsilon();
101 Eigen::Matrix < double, -1, -1, 0, MAX_MATRIX_SIZE, MAX_MATRIX_SIZE > delta =
102 (fitCov.array() - deltaCov.array()).matrix() + eps * Eigen::MatrixXd::Identity(fitCov.rows(), fitCov.cols());
103
104 fitparams.getCovariance().triangularView<Eigen::Lower>() = delta.triangularView<Eigen::Lower>();
105
106 }//end function

◆ updateState() [1/2]

void updateState ( FitParams & fitparams)

update statevector

Definition at line 66 of file KalmanCalculator.cc.

67 {
68 double eps = Eigen::NumTraits<double>::epsilon();
69 fitparams.getStateVector() = (fitparams.getStateVector().array() - (m_K * m_res).array()).matrix()
70 + eps * Eigen::MatrixXd::Identity(fitparams.getStateVector().rows(), fitparams.getStateVector().cols());
71 m_chisq = m_res.transpose() * m_Rinverse.selfadjointView<Eigen::Lower>() * m_res;
72 }

◆ updateState() [2/2]

void updateState ( FitParams & fitparams,
FitParams & oldState )

update statevector

Definition at line 74 of file KalmanCalculator.cc.

75 {
76 double eps = Eigen::NumTraits<double>::epsilon();
77 Eigen::Matrix < double, -1, 1, 0, 7, 1 > res_prime =
78 m_res + m_G * ((oldState.getStateVector().array() - fitparams.getStateVector().array()).matrix() +
79 + eps * Eigen::MatrixXd::Identity(fitparams.getStateVector().rows(), fitparams.getStateVector().cols()));
80 fitparams.getStateVector() = (oldState.getStateVector().array() - (m_K * res_prime).array()).matrix()
81 + eps * Eigen::MatrixXd::Identity(oldState.getStateVector().rows(), oldState.getStateVector().cols());
82
83 m_chisq = res_prime.transpose() * m_Rinverse.selfadjointView<Eigen::Lower>() * res_prime;
84 }

Member Data Documentation

◆ m_CGt

Eigen::Matrix< double, -1, -1, 0, MAX_MATRIX_SIZE, 7 > m_CGt
private

C times G^t.

Definition at line 88 of file KalmanCalculator.h.

◆ m_chisq

double m_chisq
private

chi2

Definition at line 59 of file KalmanCalculator.h.

◆ m_constrDim

int m_constrDim
private

dimension of the constraint

Definition at line 56 of file KalmanCalculator.h.

◆ m_G

Eigen::Matrix< double, -1, -1, 0, 7, MAX_MATRIX_SIZE > m_G
private

G former H, transforms covraince of {residuals}<->{x,p,E}.

Definition at line 76 of file KalmanCalculator.h.

◆ m_K

Eigen::Matrix< double, -1, -1, 0, MAX_MATRIX_SIZE, 7 > m_K
private

K kalman gain matrix.

Definition at line 85 of file KalmanCalculator.h.

◆ m_R

Eigen::Matrix< double, -1, -1, 0, 7, 7 > m_R
private

R residual covariance.

Definition at line 79 of file KalmanCalculator.h.

◆ m_res

Eigen::Matrix< double, -1, 1, 0, 7, 1 > m_res
private

we know the max sizes of the matrices we assume the tree is smaller than MAX_MATRIX_SIZE parameters and the largest constraint is the track constraint with 7

-> Eigen puts this on the stack

Eigen::Matrix < double, col, row, ColMajor, maxCol, maxRow> -1 = Eigen::Dynamic 0 = Eigen::ColMajor (don't touch) vector holding the residuals

Definition at line 73 of file KalmanCalculator.h.

◆ m_Rinverse

Eigen::Matrix< double, -1, -1, 0, 7, 7 > m_Rinverse
private

R inverse.

Definition at line 82 of file KalmanCalculator.h.


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