ATLAS Offline Software
AFPSiDBasicKalmanToolTrack.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2021 CERN for the benefit of the ATLAS collaboration
3 */
4 
10 
11 #ifndef AFP_LOCRECO_AFPSIDBASICKALMANTOOLTRACK_H
12 #define AFP_LOCRECO_AFPSIDBASICKALMANTOOLTRACK_H 1
13 
14 #include <list>
15 
17 
18 
19 #include <CLHEP/Matrix/Matrix.h>
20 #include <CLHEP/Matrix/Vector.h>
21 
24 {
25 public:
26 
40  AFPSiDBasicKalmanToolTrack (const xAOD::AFPSiHitsCluster* firstCluster, const xAOD::AFPSiHitsCluster* secondCluster, const CLHEP::HepMatrix& observationModel, const CLHEP::HepMatrix& observationNoise, const CLHEP::HepMatrix& covarianceMatrix);
41 
56  inline CLHEP::HepMatrix transitionModel (const double firstPos, const double secondPos) const;
57 
63  CLHEP::HepMatrix transitionModel (const double zPosition) const
64  {return transitionModel (m_zPositionHistory.back(), zPosition);}
65 
67  CLHEP::HepVector predictNextPoint (const double zPosition) const
68  {return transitionModel(zPosition)*m_positionAndSlopeHistory.back();}
69 
71  inline CLHEP::HepMatrix predictNextPointCov (const double zPosition) const
72  {CLHEP::HepMatrix transition = transitionModel(zPosition);
73  return transition*m_positionAndSlopeCovHistory.back()*transition.T();}
74 
75 
87  const xAOD::AFPSiHitsCluster* findNearestCluster (const std::vector<const xAOD::AFPSiHitsCluster*>& clusters, const float maxAllowedDistance) const;
88 
89 
91  void addCluster (const xAOD::AFPSiHitsCluster* cluster, const float htiMaxChi2);
92 
94  int holes () const {return m_holes;}
95 
97  void addHole() {++m_holes;}
98 
100  const std::list<const xAOD::AFPSiHitsCluster*>& clustersInTrack () const {return m_clustersInTrack;}
101 
103  void clearSmoothMatrices ();
104 
106  const std::list< CLHEP::HepVector >& positionAndSlopeSmooth () const {return m_positionAndSlopeSmooth;}
107 
109  const std::list< CLHEP::HepVector >& positionAndSlopeHistory () const {return m_positionAndSlopeHistory;}
110 
112  const std::list< double >& chi2Smooth() const {return m_chi2Smooth;}
113 
117  double trkChi2NDFSmooth() const;
118 
125  void smooth ();
126 protected:
127 
128 
130  std::list< float > m_zPositionHistory;
131 
137  std::list< CLHEP::HepVector > m_positionAndSlopeHistory;
138 
144  std::list< CLHEP::HepVector > m_positionAndSlopePredictedHistory;
145 
147  std::list< CLHEP::HepMatrix > m_positionAndSlopeCovHistory;
148 
150  std::list< CLHEP::HepMatrix > m_positionAndSlopeCovPredictedHistory;
151 
153  std::list< CLHEP::HepMatrix > m_transitionHistory;
154 
156  std::list< double> m_chi2History;
157 
163  std::list< CLHEP::HepVector > m_positionAndSlopeSmooth;
164 
166  std::list< CLHEP::HepMatrix > m_positionAndSlopeSmoothCov;
167 
169  std::list< double > m_chi2Smooth;
170 
172  std::list<const xAOD::AFPSiHitsCluster*> m_clustersInTrack;
173 
175  const CLHEP::HepMatrix& m_observationModel;
176 
180  const CLHEP::HepMatrix& m_observationNoise;
181 
182 
184  int m_holes;
185 };
186 
187 inline CLHEP::HepMatrix AFPSiDBasicKalmanToolTrack::transitionModel (const double firstPos, const double secondPos) const
188 {
189  CLHEP::HepMatrix transitionModel (4, 4, 1); // create unit matrix (1 on diagonal, 0 offdiagonal)
190  const double zDistance = secondPos - firstPos;
191  transitionModel[0][1] = zDistance;// set linear transformation which sets x' = x + dx*dz
192  transitionModel[2][3] = zDistance; // set linear transformation which sets y' = y + dy*dz
193 
194  return transitionModel;
195 }
196 
197 #endif
AFPSiHitsCluster.h
AFPSiDBasicKalmanToolTrack::transitionModel
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.
Definition: AFPSiDBasicKalmanToolTrack.h:187
AFPSiDBasicKalmanToolTrack::predictNextPointCov
CLHEP::HepMatrix predictNextPointCov(const double zPosition) const
Predict covariance matrix of the track at given Z.
Definition: AFPSiDBasicKalmanToolTrack.h:71
AFPSiDBasicKalmanToolTrack::m_chi2History
std::list< double > m_chi2History
History of chi2 for each step.
Definition: AFPSiDBasicKalmanToolTrack.h:156
AFPSiDBasicKalmanToolTrack::m_positionAndSlopeSmooth
std::list< CLHEP::HepVector > m_positionAndSlopeSmooth
Vectors of reconstructed true positions for each step after smoothing .
Definition: AFPSiDBasicKalmanToolTrack.h:163
AFPSiDBasicKalmanToolTrack::m_positionAndSlopeSmoothCov
std::list< CLHEP::HepMatrix > m_positionAndSlopeSmoothCov
Matrices of updated (a posteriori) estimate covariance of the position and slope for each step after ...
Definition: AFPSiDBasicKalmanToolTrack.h:166
AFPSiDBasicKalmanToolTrack::addHole
void addHole()
Increase number of holes by 1. (see m_holes)
Definition: AFPSiDBasicKalmanToolTrack.h:97
AFPSiDBasicKalmanToolTrack::m_observationModel
const CLHEP::HepMatrix & m_observationModel
Observation model, describing transformation from true to measured state.
Definition: AFPSiDBasicKalmanToolTrack.h:175
AFPSiDBasicKalmanToolTrack::m_positionAndSlopeHistory
std::list< CLHEP::HepVector > m_positionAndSlopeHistory
Vectors of reconstructed true positions for each step .
Definition: AFPSiDBasicKalmanToolTrack.h:137
AFPSiDBasicKalmanToolTrack::m_positionAndSlopeCovPredictedHistory
std::list< CLHEP::HepMatrix > m_positionAndSlopeCovPredictedHistory
Matrices of predicted estimate covariance of the position and slope for each step.
Definition: AFPSiDBasicKalmanToolTrack.h:150
AFPSiDBasicKalmanToolTrack::positionAndSlopeSmooth
const std::list< CLHEP::HepVector > & positionAndSlopeSmooth() const
Vectors of reconstructed true positions for each step after smoothing .
Definition: AFPSiDBasicKalmanToolTrack.h:106
AFPSiDBasicKalmanToolTrack::m_clustersInTrack
std::list< const xAOD::AFPSiHitsCluster * > m_clustersInTrack
Clusters used to reconstruct the track.
Definition: AFPSiDBasicKalmanToolTrack.h:172
AFPSiDBasicKalmanToolTrack::clearSmoothMatrices
void clearSmoothMatrices()
Clears matrices used for smoothing the track.
Definition: AFPSiDBasicKalmanToolTrack.cxx:115
AFPSiDBasicKalmanToolTrack::trkChi2NDFSmooth
double trkChi2NDFSmooth() const
Chi2 per degrees of freedom (clusters) for the whole smoothed track.
Definition: AFPSiDBasicKalmanToolTrack.cxx:216
AFPSiDBasicKalmanToolTrack::m_zPositionHistory
std::list< float > m_zPositionHistory
Vectors of reconstructed true positions for each step.
Definition: AFPSiDBasicKalmanToolTrack.h:130
AFPSiDBasicKalmanToolTrack::clustersInTrack
const std::list< const xAOD::AFPSiHitsCluster * > & clustersInTrack() const
Clusters used to reconstruct the track.
Definition: AFPSiDBasicKalmanToolTrack.h:100
AFPSiDBasicKalmanToolTrack::findNearestCluster
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.
Definition: AFPSiDBasicKalmanToolTrack.cxx:97
AFPSiDBasicKalmanToolTrack::addCluster
void addCluster(const xAOD::AFPSiHitsCluster *cluster, const float htiMaxChi2)
Adds a new cluster to the track and updates history and position.
Definition: AFPSiDBasicKalmanToolTrack.cxx:48
AFPSiDBasicKalmanToolTrack::m_positionAndSlopeCovHistory
std::list< CLHEP::HepMatrix > m_positionAndSlopeCovHistory
Matrices of updated (a posteriori) estimate covariance of the position and slope for each step .
Definition: AFPSiDBasicKalmanToolTrack.h:147
AFPSiDBasicKalmanToolTrack::transitionModel
CLHEP::HepMatrix transitionModel(const double zPosition) const
Create state transition matrix to calculate x and y position and slopes at given Z.
Definition: AFPSiDBasicKalmanToolTrack.h:63
AFPSiDBasicKalmanToolTrack::m_positionAndSlopePredictedHistory
std::list< CLHEP::HepVector > m_positionAndSlopePredictedHistory
Vectors of predicted true position of a track for each step.
Definition: AFPSiDBasicKalmanToolTrack.h:144
AFPSiDBasicKalmanToolTrack
Class representing a reconstructed track with Kalman filter.
Definition: AFPSiDBasicKalmanToolTrack.h:24
AFPSiDBasicKalmanToolTrack::m_transitionHistory
std::list< CLHEP::HepMatrix > m_transitionHistory
Matrices of transitionModel used at each step of track reconstruction.
Definition: AFPSiDBasicKalmanToolTrack.h:153
AFPSiDBasicKalmanToolTrack::holes
int holes() const
Definition: AFPSiDBasicKalmanToolTrack.h:94
AFPSiDBasicKalmanToolTrack::positionAndSlopeHistory
const std::list< CLHEP::HepVector > & positionAndSlopeHistory() const
Vectors of reconstructed true positions for each step .
Definition: AFPSiDBasicKalmanToolTrack.h:109
AFPSiDBasicKalmanToolTrack::m_holes
int m_holes
Number of layers in which hit from a track is expected, but is not present.
Definition: AFPSiDBasicKalmanToolTrack.h:184
AFPSiDBasicKalmanToolTrack::predictNextPoint
CLHEP::HepVector predictNextPoint(const double zPosition) const
Predict x and y positions and slopes of the track at given Z.
Definition: AFPSiDBasicKalmanToolTrack.h:67
AFPSiDBasicKalmanToolTrack::smooth
void smooth()
Run smoothing algorithm.
Definition: AFPSiDBasicKalmanToolTrack.cxx:122
RunTileMonitoring.clusters
clusters
Definition: RunTileMonitoring.py:133
AFPSiDBasicKalmanToolTrack::chi2Smooth
const std::list< double > & chi2Smooth() const
History of chi2 for each step after smoothing.
Definition: AFPSiDBasicKalmanToolTrack.h:112
xAOD::AFPSiHitsCluster_v1
Class representing a cluster of AFP pixel hits.
Definition: AFPSiHitsCluster_v1.h:32
AFPSiDBasicKalmanToolTrack::m_observationNoise
const CLHEP::HepMatrix & m_observationNoise
Matrix representing observation noise disturbing the measurement.
Definition: AFPSiDBasicKalmanToolTrack.h:180
AFPSiDBasicKalmanToolTrack::m_chi2Smooth
std::list< double > m_chi2Smooth
History of chi2 for each step after smoothing.
Definition: AFPSiDBasicKalmanToolTrack.h:169
AFPSiDBasicKalmanToolTrack::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.
Definition: AFPSiDBasicKalmanToolTrack.cxx:19