ATLAS Offline Software
|
Class representing a reconstructed track with Kalman filter. More...
#include <AFPSiDBasicKalmanToolTrack.h>
Public Member Functions | |
AFPSiDBasicKalmanToolTrack (const xAOD::AFPSiHitsCluster *firstCluster, const xAOD::AFPSiHitsCluster *secondCluster, const CLHEP::HepMatrix &observationModel, const CLHEP::HepMatrix &observationNoise, const CLHEP::HepMatrix &covarianceMatrix) | |
Main constructor that creates the track seed from two clusters. More... | |
CLHEP::HepMatrix | transitionModel (const double firstPos, const double secondPos) const |
Create state transition matrix to calculate x and y position and slopes between to Z positions. More... | |
CLHEP::HepMatrix | transitionModel (const double zPosition) const |
Create state transition matrix to calculate x and y position and slopes at given Z. More... | |
CLHEP::HepVector | predictNextPoint (const double zPosition) const |
Predict x and y positions and slopes of the track at given Z. More... | |
CLHEP::HepMatrix | predictNextPointCov (const double zPosition) const |
Predict covariance matrix of the track at given Z. More... | |
const xAOD::AFPSiHitsCluster * | findNearestCluster (const std::vector< const xAOD::AFPSiHitsCluster * > &clusters, const float maxAllowedDistance) const |
Finds the cluster which is nearest to the track, returns nullptr if no is found. More... | |
void | addCluster (const xAOD::AFPSiHitsCluster *cluster, const float htiMaxChi2) |
Adds a new cluster to the track and updates history and position. More... | |
int | holes () const |
void | addHole () |
Increase number of holes by 1. (see m_holes) More... | |
const std::list< const xAOD::AFPSiHitsCluster * > & | clustersInTrack () const |
Clusters used to reconstruct the track. More... | |
void | clearSmoothMatrices () |
Clears matrices used for smoothing the track. More... | |
const std::list< CLHEP::HepVector > & | positionAndSlopeSmooth () const |
Vectors of reconstructed true positions for each step after smoothing \(x_{k}\). More... | |
const std::list< CLHEP::HepVector > & | positionAndSlopeHistory () const |
Vectors of reconstructed true positions for each step \(x_{k}\). More... | |
const std::list< double > & | chi2Smooth () const |
History of chi2 for each step after smoothing. More... | |
double | trkChi2NDFSmooth () const |
Chi2 per degrees of freedom (clusters) for the whole smoothed track. More... | |
void | smooth () |
Run smoothing algorithm. More... | |
Protected Attributes | |
std::list< float > | m_zPositionHistory |
Vectors of reconstructed true positions for each step. More... | |
std::list< CLHEP::HepVector > | m_positionAndSlopeHistory |
Vectors of reconstructed true positions for each step \(x_{k}\). More... | |
std::list< CLHEP::HepVector > | m_positionAndSlopePredictedHistory |
Vectors of predicted true position of a track for each step. More... | |
std::list< CLHEP::HepMatrix > | m_positionAndSlopeCovHistory |
Matrices of updated (a posteriori) estimate covariance of the position and slope for each step \(P_{k|k}\). More... | |
std::list< CLHEP::HepMatrix > | m_positionAndSlopeCovPredictedHistory |
Matrices of predicted estimate covariance of the position and slope for each step. More... | |
std::list< CLHEP::HepMatrix > | m_transitionHistory |
Matrices of transitionModel used at each step of track reconstruction. More... | |
std::list< double > | m_chi2History |
History of chi2 for each step. More... | |
std::list< CLHEP::HepVector > | m_positionAndSlopeSmooth |
Vectors of reconstructed true positions for each step after smoothing \(x_{k}\). More... | |
std::list< CLHEP::HepMatrix > | m_positionAndSlopeSmoothCov |
Matrices of updated (a posteriori) estimate covariance of the position and slope for each step after smoothing \(P_{k|k}\). More... | |
std::list< double > | m_chi2Smooth |
History of chi2 for each step after smoothing. More... | |
std::list< const xAOD::AFPSiHitsCluster * > | m_clustersInTrack |
Clusters used to reconstruct the track. More... | |
const CLHEP::HepMatrix & | m_observationModel |
Observation model, describing transformation from true to measured state. More... | |
const CLHEP::HepMatrix & | m_observationNoise |
Matrix representing observation noise disturbing the measurement. More... | |
int | m_holes |
Number of layers in which hit from a track is expected, but is not present. More... | |
Class representing a reconstructed track with Kalman filter.
Definition at line 23 of file AFPSiDBasicKalmanToolTrack.h.
AFPSiDBasicKalmanToolTrack::AFPSiDBasicKalmanToolTrack | ( | const xAOD::AFPSiHitsCluster * | firstCluster, |
const xAOD::AFPSiHitsCluster * | secondCluster, | ||
const CLHEP::HepMatrix & | observationModel, | ||
const CLHEP::HepMatrix & | observationNoise, | ||
const CLHEP::HepMatrix & | covarianceMatrix | ||
) |
Main constructor that creates the track seed from two clusters.
It adds the two clusters to the list of clusters used to reconstruct the track. It sets the position of the last point of the track to be equal to the position of the second cluster and calculates track slopes using first and second clusters. It adds first entry to m_chi2History and sets it to 0.
firstCluster | pointer to the first cluster of the track |
secondCluster | pointer to the second cluster of the track |
observationModel | matrix representing observation model (see m_observationModel) |
observationNoise | see m_observationNoise |
covarianceMatrix | initial error covariance matrix (see m_positionAndSlopeCovHistory) |
Definition at line 19 of file AFPSiDBasicKalmanToolTrack.cxx.
void AFPSiDBasicKalmanToolTrack::addCluster | ( | const xAOD::AFPSiHitsCluster * | cluster, |
const float | htiMaxChi2 | ||
) |
Adds a new cluster to the track and updates history and position.
Definition at line 48 of file AFPSiDBasicKalmanToolTrack.cxx.
|
inline |
Increase number of holes by 1. (see m_holes)
Definition at line 97 of file AFPSiDBasicKalmanToolTrack.h.
|
inline |
History of chi2 for each step after smoothing.
Definition at line 112 of file AFPSiDBasicKalmanToolTrack.h.
void AFPSiDBasicKalmanToolTrack::clearSmoothMatrices | ( | ) |
Clears matrices used for smoothing the track.
Definition at line 115 of file AFPSiDBasicKalmanToolTrack.cxx.
|
inline |
Clusters used to reconstruct the track.
Definition at line 100 of file AFPSiDBasicKalmanToolTrack.h.
const xAOD::AFPSiHitsCluster * AFPSiDBasicKalmanToolTrack::findNearestCluster | ( | const std::vector< const xAOD::AFPSiHitsCluster * > & | clusters, |
const float | maxAllowedDistance | ||
) | const |
Finds the cluster which is nearest to the track, returns nullptr if no is found.
The method predicts what are track X and Y coordinates for the Z of a cluster. Next, loops over the clusters in the vector and calculates distance \(r = \sqrt{\Delta x^2 + \Delta y^2}\). Returns nullptr if there is no cluster which is closer than #maxAllowedDistance or if there are clusters with r < #maxAllowedDistance a pointer to the closest one.
clusters | vector of clusters from which the closest one is to be found |
maxAllowedDistance | if no cluster is closer than this argument nullptr is returned |
Definition at line 97 of file AFPSiDBasicKalmanToolTrack.cxx.
|
inline |
Definition at line 94 of file AFPSiDBasicKalmanToolTrack.h.
|
inline |
Vectors of reconstructed true positions for each step \(x_{k}\).
Each 4D vector containing x, y position of the track in a given step as well as slopes in these two directions \(\Delta x, \Delta y\). \(x_{k} = (x, \Delta x, y, \Delta y)\)
Definition at line 109 of file AFPSiDBasicKalmanToolTrack.h.
|
inline |
Vectors of reconstructed true positions for each step after smoothing \(x_{k}\).
Each 4D vector containing x, y position of the track in a given step as well as slopes in these two directions \(\Delta x, \Delta y\). \(x_{k} = (x, \Delta x, y, \Delta y)\)
Definition at line 106 of file AFPSiDBasicKalmanToolTrack.h.
|
inline |
Predict x and y positions and slopes of the track at given Z.
Definition at line 67 of file AFPSiDBasicKalmanToolTrack.h.
|
inline |
Predict covariance matrix of the track at given Z.
Definition at line 71 of file AFPSiDBasicKalmanToolTrack.h.
void AFPSiDBasicKalmanToolTrack::smooth | ( | ) |
Run smoothing algorithm.
The algorithm should probably just reconstruct track starting from the last hit, but the covaraince matrices are changed, so the result maybe slightly different. Algorithm is just copied from the previous implementation.
Definition at line 122 of file AFPSiDBasicKalmanToolTrack.cxx.
|
inline |
Create state transition matrix to calculate x and y position and slopes between to Z positions.
The matrix leaves slopes unchanged and adds to the previous positions slopes multiplied by the z distance i.e. secondPos - firstPos.
\[ \left( \begin{array}{cccc} 1 & \Delta z & 0 & 0 0 & 1 & 0 & 0 0 & 0 & 1 & \Delta z 0 & 0 & 0 & 1 \end{array} \right) \]
Definition at line 187 of file AFPSiDBasicKalmanToolTrack.h.
|
inline |
Create state transition matrix to calculate x and y position and slopes at given Z.
The matrix for transition between the last point and given Z position is calculated.
Definition at line 63 of file AFPSiDBasicKalmanToolTrack.h.
double AFPSiDBasicKalmanToolTrack::trkChi2NDFSmooth | ( | ) | const |
Chi2 per degrees of freedom (clusters) for the whole smoothed track.
Sum chi2 for each cluster and divide by the number of clusters.
Definition at line 216 of file AFPSiDBasicKalmanToolTrack.cxx.
|
protected |
History of chi2 for each step.
Definition at line 156 of file AFPSiDBasicKalmanToolTrack.h.
|
protected |
History of chi2 for each step after smoothing.
Definition at line 169 of file AFPSiDBasicKalmanToolTrack.h.
|
protected |
Clusters used to reconstruct the track.
Definition at line 172 of file AFPSiDBasicKalmanToolTrack.h.
|
protected |
Number of layers in which hit from a track is expected, but is not present.
Definition at line 184 of file AFPSiDBasicKalmanToolTrack.h.
|
protected |
Observation model, describing transformation from true to measured state.
Definition at line 175 of file AFPSiDBasicKalmanToolTrack.h.
|
protected |
Matrix representing observation noise disturbing the measurement.
Used to estimate covariance matrix of the track.
Definition at line 180 of file AFPSiDBasicKalmanToolTrack.h.
|
protected |
Matrices of updated (a posteriori) estimate covariance of the position and slope for each step \(P_{k|k}\).
Definition at line 147 of file AFPSiDBasicKalmanToolTrack.h.
|
protected |
Matrices of predicted estimate covariance of the position and slope for each step.
Definition at line 150 of file AFPSiDBasicKalmanToolTrack.h.
|
protected |
Vectors of reconstructed true positions for each step \(x_{k}\).
Each 4D vector containing x, y position of the track in a given step as well as slopes in these two directions \(\Delta x, \Delta y\). \(x_{k} = (x, \Delta x, y, \Delta y)\)
Definition at line 137 of file AFPSiDBasicKalmanToolTrack.h.
|
protected |
Vectors of predicted true position of a track for each step.
Each 4D vector containing x, y position of the track in a given step as well as slopes in these two directions \(\Delta x, \Delta y\). \(x_{k} = (x, \Delta x, y, \Delta y)\)
Definition at line 144 of file AFPSiDBasicKalmanToolTrack.h.
|
protected |
Vectors of reconstructed true positions for each step after smoothing \(x_{k}\).
Each 4D vector containing x, y position of the track in a given step as well as slopes in these two directions \(\Delta x, \Delta y\). \(x_{k} = (x, \Delta x, y, \Delta y)\)
Definition at line 163 of file AFPSiDBasicKalmanToolTrack.h.
|
protected |
Matrices of updated (a posteriori) estimate covariance of the position and slope for each step after smoothing \(P_{k|k}\).
Definition at line 166 of file AFPSiDBasicKalmanToolTrack.h.
|
protected |
Matrices of transitionModel used at each step of track reconstruction.
Definition at line 153 of file AFPSiDBasicKalmanToolTrack.h.
|
protected |
Vectors of reconstructed true positions for each step.
Definition at line 130 of file AFPSiDBasicKalmanToolTrack.h.