ATLAS Offline Software
AFPSiDBasicKalmanToolTrack.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3 */
4 
10 
11 #include <cmath>
12 
13 // debug remove after debugging
14 #include<iostream>
15 // end remove after debugging
16 
18 
19 AFPSiDBasicKalmanToolTrack::AFPSiDBasicKalmanToolTrack (const xAOD::AFPSiHitsCluster* firstCluster, const xAOD::AFPSiHitsCluster* secondCluster, const CLHEP::HepMatrix& observationModel, const CLHEP::HepMatrix& observationNoise, const CLHEP::HepMatrix& covarianceMatrix) :
20  m_observationModel(observationModel),
21  m_observationNoise(observationNoise),
22  m_holes(0)
23 {
24  m_positionAndSlopeHistory.emplace_back(4);
25 
26  const double zDistance = secondCluster->zLocal() - firstCluster->zLocal(); // absolute value taken in case global z position is used, in that case the distance would be negative for stations on the A side
27 
28  CLHEP::HepVector& positionAndSlope = m_positionAndSlopeHistory.back();
29  positionAndSlope[0] = secondCluster->xLocal(); // set track X position
30  positionAndSlope[1] = (secondCluster->xLocal() - firstCluster->xLocal())/zDistance; // set track X slope
31  positionAndSlope[2] = secondCluster->yLocal(); // set track Y position
32  positionAndSlope[3] = (secondCluster->yLocal() - firstCluster->yLocal())/zDistance; // set track Y slope
33 
34  // set z position of the track
35  m_zPositionHistory.push_back(secondCluster->zLocal());
36 
37  // add clusters to the list
38  m_clustersInTrack.push_back(firstCluster);
39  m_clustersInTrack.push_back(secondCluster);
40 
41  // set the first covariance matrix
42  m_positionAndSlopeCovHistory.emplace_back(covarianceMatrix);
43 
44  // set chi2 to 0 for two points
45  m_chi2History.emplace_back(0);
46 }
47 
48 void AFPSiDBasicKalmanToolTrack::addCluster (const xAOD::AFPSiHitsCluster* cluster, const float clusterMaxChi2)
49 {
50  // prepare a vector with new observation point coordinates
51  CLHEP::HepVector newPoint (2);
52  newPoint[0] = cluster->xLocal();
53  newPoint[1] = cluster->yLocal();
54 
55  // get the values for the predicted new point
56  const double newZ = cluster->zLocal();
57 
58  const CLHEP::HepVector predictedTruePosition = predictNextPoint(newZ);
59  const CLHEP::HepMatrix predictedCovariance = predictNextPointCov(newZ);
60 
61  // difference between measured and predicted (called innovation y_k and S_k)
62  const CLHEP::HepVector residualPredicted = newPoint - m_observationModel*predictedTruePosition;
63  const CLHEP::HepMatrix residualPredictedCov = m_observationNoise + m_observationModel*predictedCovariance*m_observationModel.T();
64 
65  // optimal Kalman gain (Kk)
66  const CLHEP::HepMatrix kalmanGain = predictedCovariance*m_observationModel.T()*qr_inverse(residualPredictedCov);
67 
68  // === calculated new values ===
69  // Updated (a posteriori) state estimate
70  const CLHEP::HepVector newPosition = predictedTruePosition + kalmanGain*residualPredicted;
71 
72  // Updated (a posteriori) residuals
73  const CLHEP::HepVector newResidual = newPoint - m_observationModel*newPosition;
74  const CLHEP::HepMatrix newResidualCov = m_observationNoise - m_observationModel*kalmanGain*m_observationNoise;
75 
76  // chi2 statistics
77  const CLHEP::HepVector chi2Vector = newResidual.T()*qr_inverse(newResidualCov)*newResidual;
78 
79  // check if hit is good
80  if (chi2Vector[0] < clusterMaxChi2) {
81  // add cluster and all information
82  m_clustersInTrack.push_back(cluster);
83  m_zPositionHistory.push_back(newZ);
84  m_transitionHistory.emplace_back(transitionModel(newZ));
85  m_positionAndSlopePredictedHistory.push_back(predictedTruePosition);
86  m_positionAndSlopeCovPredictedHistory.push_back(predictedCovariance);
87  m_positionAndSlopeHistory.push_back(newPosition);
88  m_chi2History.emplace_back(chi2Vector[0]);
89  // Updated (a posteriori) estimate covariance
90  m_positionAndSlopeCovHistory.emplace_back(predictedCovariance - kalmanGain*m_observationModel*predictedCovariance);
91  }
92  else {
93  addHole();
94  }
95 }
96 
97 const xAOD::AFPSiHitsCluster* AFPSiDBasicKalmanToolTrack::findNearestCluster (const std::vector<const xAOD::AFPSiHitsCluster*>& clusters, const float maxAllowedDistance) const
98 {
99  double minDisntace = maxAllowedDistance;
100  const xAOD::AFPSiHitsCluster* nearestCluster = nullptr;
101  for (const xAOD::AFPSiHitsCluster* theCluster : clusters) {
102  const CLHEP::HepVector predictedPoint = predictNextPoint(theCluster->zLocal());
103  const double xDiff = predictedPoint[0] - theCluster->xLocal();
104  const double yDiff = predictedPoint[2] - theCluster->yLocal(); // index 2, because this is (x, dx, y, dy) vector
105  const double distance = std::sqrt(xDiff*xDiff + yDiff*yDiff);
106  if (distance < minDisntace) {
107  minDisntace = distance;
108  nearestCluster = theCluster;
109  }
110  }
111 
112  return nearestCluster;
113 }
114 
116 {
117  m_positionAndSlopeSmooth.clear();
119  m_chi2Smooth.clear();
120 }
121 
123 {
125 
126  // prepare iterators for going backwards through the lists
127  std::list< CLHEP::HepVector >::const_iterator positionIter = m_positionAndSlopeHistory.end();
128  std::list< CLHEP::HepVector >::const_iterator positionBegin = m_positionAndSlopeHistory.begin();
129  // if no position is reconstructed there is nothing to be done
130  if (positionIter != positionBegin) {
131  // prepare remaining iterators
132  std::list< CLHEP::HepMatrix >::const_iterator covIter = m_positionAndSlopeCovHistory.end();
133  std::list< const xAOD::AFPSiHitsCluster* >::const_iterator clustersIter = m_clustersInTrack.end();
134  std::list< CLHEP::HepVector >::const_iterator positionPredictIter = m_positionAndSlopePredictedHistory.end();
135  std::list< CLHEP::HepMatrix >::const_iterator covPredictIter = m_positionAndSlopeCovPredictedHistory.end();
136  std::list< CLHEP::HepMatrix >::const_iterator transitionIter = m_transitionHistory.end();
137 
138  // === prepare first point ===
139  // iterator moved back by 1, because end() points after the last element
140  const CLHEP::HepVector& lastPosition = *(--positionIter);
141  const CLHEP::HepMatrix& lastCovariance = *(--covIter);
142  const xAOD::AFPSiHitsCluster* lastCluster = *(--clustersIter);
143 
144  m_positionAndSlopeSmooth.push_back(lastPosition);
145  m_positionAndSlopeSmoothCov.push_back(lastCovariance);
146 
147  // prepare a vector with new observation point coordinates
148  CLHEP::HepVector lastPoint (2);
149  lastPoint[0] = lastCluster->xLocal();
150  lastPoint[1] = lastCluster->yLocal();
151 
152  const CLHEP::HepVector residualSmooth = lastPoint - m_observationModel*lastPosition;
153  const CLHEP::HepMatrix residualSmoothCov = m_observationNoise - m_observationModel*lastCovariance*m_observationModel.T();
154  const CLHEP::HepVector chi2Smooth = residualSmooth.T()*qr_inverse(residualSmoothCov)*residualSmooth;
155  m_chi2Smooth.push_back(chi2Smooth[0]);
156 
157  // === loop over remaining clusters ===
158  while (positionIter != positionBegin) {
159  // --- read needed variables ---
160  const CLHEP::HepVector& position = *(--positionIter);
161  const CLHEP::HepMatrix& covariance = *(--covIter);
162  const xAOD::AFPSiHitsCluster* cluster = *(--clustersIter);
163 
164  // shifted by 1 in history with respect to the position and covariance
165  const CLHEP::HepVector& prevPositionPred = *(--positionPredictIter);
166  const CLHEP::HepMatrix& prevCovPred = *(--covPredictIter);
167  const CLHEP::HepMatrix& prevTransition = *(--transitionIter);
168 
169  // --- calculate the smoothed position and covariance ---
170  const CLHEP::HepMatrix Ai = qr_inverse(prevTransition); // when Qk == 0
171  // const CLHEP::HepMatrix Ai = covariance * prevTransition.T() * qr_inverse(prevCovPred);
172  const CLHEP::HepVector positionSmooth = position + Ai*( m_positionAndSlopeSmooth.back() - prevPositionPred );
173  const CLHEP::HepMatrix covarianceSmooth = covariance + Ai*( m_positionAndSlopeSmoothCov.back() - prevCovPred )*Ai.T();
174 
175  // --- calculate smoothed chi2 ---
176  CLHEP::HepVector point (2);
177  point[0] = cluster->xLocal();
178  point[1] = cluster->yLocal();
179  const CLHEP::HepVector residualSmooth = point - m_observationModel*positionSmooth;
180  const CLHEP::HepMatrix residualCovSmooth = m_observationNoise - m_observationModel*covarianceSmooth*m_observationModel.T();
181  CLHEP::HepVector chi2Smooth = residualSmooth.T()*qr_inverse(residualCovSmooth)*residualSmooth;
182 
183  // --- add results to the lists ---
184  m_positionAndSlopeSmooth.push_back(positionSmooth);
185  m_positionAndSlopeSmoothCov.push_back(covarianceSmooth);
186  m_chi2Smooth.push_back(chi2Smooth[0]);
187  } // end while over points/clusters
188 
189 
190  // === extrapolate to the first plane ===
191  clustersIter = m_clustersInTrack.begin();
192  const xAOD::AFPSiHitsCluster* firstCluster = *clustersIter;
193  const xAOD::AFPSiHitsCluster* secondCluster = *(++clustersIter);
194  const CLHEP::HepMatrix firstTransition = transitionModel(firstCluster->zLocal(), secondCluster->zLocal());
195 
196  const CLHEP::HepMatrix Ai = qr_inverse(firstTransition); // when Qk == 0
197  const CLHEP::HepVector firstPosition = Ai*(m_positionAndSlopeSmooth.back());
198  const CLHEP::HepMatrix firstCovariance = m_positionAndSlopeCovHistory.front() + Ai*( m_positionAndSlopeSmoothCov.back() - m_positionAndSlopeCovHistory.front())*Ai.T();
199 
200  // --- calculate smoothed chi2 ---
201  const CLHEP::HepVector firstResidual = m_observationModel*Ai*(m_positionAndSlopeHistory.front() - m_positionAndSlopeSmooth.back());
202  CLHEP::HepVector firstPoint (2);
203  firstPoint[0] = firstCluster->xLocal();
204  firstPoint[1] = firstCluster->yLocal();
205  // const CLHEP::HepVector firstResidual = firstPoint - m_observationModel*firstPosition;
206  const CLHEP::HepMatrix firstResidualCov = m_observationNoise - m_observationModel*firstCovariance*m_observationModel.T();
207  CLHEP::HepVector firstChi2 = firstResidual.T()*qr_inverse(firstResidualCov)*firstResidual;
208 
209  // --- add results to the lists ---
210  m_positionAndSlopeSmooth.push_back(firstPosition);
211  m_positionAndSlopeSmoothCov.push_back(firstCovariance);
212  m_chi2Smooth.push_back(firstChi2[0]);
213  }
214 }
215 
217 {
218  double chi2 = 0;
219  for (double val : m_chi2Smooth) {
220  chi2 += val;
221  }
222 
223  return chi2/m_chi2Smooth.size();
224 }
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::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.h
Header file for AFPSiDBasicKalmanToolTrack used in tracks reconstruction with Kalman filter.
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
chi2
double chi2(TH1 *h0, TH1 *h1)
Definition: comparitor.cxx:522
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
xAOD::AFPSiHitsCluster_v1::xLocal
float xLocal() const
Cluster position along X axis in station local coordinate system.
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::m_positionAndSlopePredictedHistory
std::list< CLHEP::HepVector > m_positionAndSlopePredictedHistory
Vectors of predicted true position of a track for each step.
Definition: AFPSiDBasicKalmanToolTrack.h:144
AFPSiDBasicKalmanToolTrack::m_transitionHistory
std::list< CLHEP::HepMatrix > m_transitionHistory
Matrices of transitionModel used at each step of track reconstruction.
Definition: AFPSiDBasicKalmanToolTrack.h:153
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
Pythia8_RapidityOrderMPI.val
val
Definition: Pythia8_RapidityOrderMPI.py:14
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
xAOD::AFPSiHitsCluster_v1::yLocal
float yLocal() const
Cluster position along Y axis in station local coordinate system.
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
Amg::distance
float distance(const Amg::Vector3D &p1, const Amg::Vector3D &p2)
calculates the distance between two point in 3D space
Definition: GeoPrimitivesHelpers.h:54
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
xAOD::AFPSiHitsCluster_v1::zLocal
float zLocal() const
Cluster position along Z axis in station local coordinate system.