ATLAS Offline Software
GsfMeasurementUpdator.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3 */
11 // Suppress a bogus warning.
12 #if __GNUC__ >= 14 && defined(__aarch64__)
13 # pragma GCC diagnostic ignored "-Warray-bounds"
14 #endif
15 
17 //
21 //
25 
26 #include "CxxUtils/inline_hints.h"
27 
28 #include <memory>
29 
30 namespace {
31 
32 // constants
33 const Trk::ProjectionMatricesSet s_reMatrices(5);
34 constexpr double s_thetaGainDampingValue = 0.1;
35 const AmgVector(5) s_cov0Vec {250., 250., 0.25, 0.25, 0.000001};
36 const AmgMatrix(5, 5) s_unitMatrix(AmgMatrix(5, 5)::Identity());
37 const AmgVector(2) s_thetaMin(0.0, -M_PI);
38 enum RangeCheckDef
39 {
40  absoluteCheck = 0,
42 };
43 
44 struct componentsCache
45 {
46  struct element
47  {
48  double determinantR;
49  double chi2;
50  };
51  std::array<element, GSFConstants::maxNumberofStateComponents> elements{};
52  size_t numElements = 0;
53 };
54 
56 int measurementCoord_1D(int paramKey) {
57  for (int i = 0; i < 5; ++i) {
58  if (paramKey & (1 << i)) {
59  return i;
60  }
61  }
62 
63  return 0;
64 }
65 
72 inline bool
73 thetaPhiWithinRange_5D(const AmgVector(5) & V, const RangeCheckDef rcd)
74 {
75  return ((std::abs(V(Trk::phi)) <= M_PI) &&
76  (V(Trk::theta) >= s_thetaMin(rcd)) && (V(Trk::theta) <= M_PI));
77 }
78 
82 inline bool
83 thetaWithinRange_5D(const AmgVector(5) & V)
84 {
85  return (V(Trk::theta) >= 0.0 && (V(Trk::theta) <= M_PI));
86 }
87 
88 bool
89 correctThetaPhiRange_5D(AmgVector(5) & V,
90  AmgSymMatrix(5) & M,
91  const RangeCheckDef rcd)
92 {
93 
94  // correct theta coordinate
95  if (V(Trk::theta) < s_thetaMin((int)rcd) || V(Trk::theta) > M_PI) {
96  // absolute theta: repair if between -pi and +2pi.
97  // differential theta: repair if between -pi and +pi
98  if (V(Trk::theta) < -M_PI) {
99  return false;
100  }
101  if (V(Trk::theta) > (rcd == differentialCheck ? M_PI : 2 * M_PI)) {
102  return false;
103  }
104  if (V(Trk::theta) > M_PI) {
105  V(Trk::theta) = 2 * M_PI - V(Trk::theta);
106  V(Trk::phi) += (V(Trk::phi) > 0.0) ? -M_PI : M_PI;
107  } else if (V(Trk::theta) < 0.0) {
108  V(Trk::theta) = -V(Trk::theta);
109  V(Trk::phi) += (V(Trk::phi) > 0.0) ? -M_PI : M_PI;
110  }
111  // correct also cov matrix:
112  M.fillSymmetric(0, 3, -M(0, 3)); // cov(polar,locX)
113  M.fillSymmetric(1, 3, -M(1, 3)); // cov(polar,locY)
114  M.fillSymmetric(2, 3, -M(2, 3)); // cov(polar,azimuthal)
115  M.fillSymmetric(3, 4, -M(3, 4)); // cov(polar,qOverP)
116  }
117 
118  // correct phi coordinate if necessary
119  if ((V(Trk::phi) > M_PI)) {
120  V(Trk::phi) = fmod(V(Trk::phi) + M_PI, 2 * M_PI) - M_PI;
121  } else if ((V(Trk::phi) < -M_PI)) {
122  V(Trk::phi) = fmod(V(Trk::phi) - M_PI, 2 * M_PI) + M_PI;
123  }
124 
125  return true;
126 }
127 
132 template<int DIM>
133 AmgSymMatrix(DIM) projection_T(const AmgSymMatrix(5) & M, int key)
134 {
135  if (key == 3 || key == 7 || key == 15) { // shortcuts for the most common use
136  // cases
137  return M.topLeftCorner<DIM, DIM>();
138  }
139 
140  Eigen::Matrix<int, DIM, 1, 0, DIM, 1> iv;
141  iv.setZero();
142  for (int i = 0, k = 0; i < 5; ++i) {
143  iv[k++] = key & (1 << i) ? i : 0;
144  }
145 
146  AmgSymMatrix(DIM) covSubMatrix;
147  covSubMatrix.setZero();
148  for (int i = 0; i < DIM; ++i) {
149  for (int j = 0; j < DIM; ++j) {
150  covSubMatrix(i, j) = M(iv(i), iv(j));
151  }
152  }
153 
154  return covSubMatrix;
155 }
156 
165 bool
166 calculateFilterStep_1D(Trk::TrackParameters& TP,
167  const AmgSymMatrix(5) & trkCov,
168  double measPar,
169  double measCov,
170  int paramKey,
171  Trk::FitQualityOnSurface& fitQoS)
172 {
173  const int mk = measurementCoord_1D(paramKey);
174  // get the parameters from the
175  const AmgVector(5)& trkPar = TP.parameters();
176  // use measuring coordinate (variable "mk") instead of reduction matrix
177  const double r = measPar - trkPar(mk);
178  double R = measCov + trkCov(mk, mk);
179  if (R == 0.0) {
180  return false;
181  }
182  R = 1. / R;
183  // --- compute Kalman gain matrix
184  AmgVector(5) K = trkCov.col(mk) * R;
185  // --- compute local filtered state, here = TP+K*r = TP + TCov * H.T * R * r
186  AmgVector(5) newPar = trkPar + K * r;
187 
188  if (!thetaWithinRange_5D(newPar)) {
189 
190  if (mk != Trk::theta &&
191  (std::abs(R * r) > 1.0 ||
192  trkCov(Trk::theta, Trk::theta) > 0.1 * s_cov0Vec(Trk::theta))) {
193 
194  AmgVector(5) dampedCov = trkCov.col(mk);
195 
196  dampedCov(Trk::theta) = dampedCov(Trk::theta) * s_thetaGainDampingValue;
197 
198  newPar = trkPar + dampedCov * R * r;
199 
200  K(Trk::theta, 0) = K(Trk::theta, 0) * s_thetaGainDampingValue;
201  }
202  }
203  // --- compute covariance matrix of local filteres state
204  AmgMatrix(5, 5) KtimesH;
205  KtimesH.setZero();
206  KtimesH.col(mk) = K;
207  const AmgMatrix(5, 5) M = s_unitMatrix - KtimesH;
208  AmgSymMatrix(5) newCov =
209  trkCov.similarity(M) + K * measCov * K.transpose();
210 
211  if (!thetaPhiWithinRange_5D(newPar, absoluteCheck) &&
212  !correctThetaPhiRange_5D(newPar, newCov, absoluteCheck)) {
213  return false;
214  }
215 
216  const double predictedResidual = (measPar - newPar(mk));
217  const AmgSymMatrix(5)& updatedCov = newCov;
218 
219  double chiSquared = measCov - updatedCov(mk, mk);
220  if (chiSquared != 0.0) {
221  // get chi2 = r.T() * R^-1 * r
222  chiSquared = predictedResidual * predictedResidual / chiSquared;
223  }
225  // In place update of parameters
226  TP.updateParameters(newPar, newCov);
227  return true;
228 }
229 
230 bool
231 calculateFilterStep_5D(Trk::TrackParameters& TP,
232  const AmgSymMatrix(5) & trkCov,
233  const AmgVector(5) & measPar,
234  const AmgSymMatrix(5) & measCov,
236 {
237  // get the parameter vector
238  const AmgVector(5)& trkPar = TP.parameters();
239  // Kalman gain K, residual r, combined covariance R
240  const AmgVector(5) r = measPar - trkPar;
241  // for full safety in Eigen see
242  // http://eigen.tuxfamily.org/dox/classEigen_1_1FullPivLU.html
243  const AmgSymMatrix(5) R = (measCov + trkCov).inverse();
244  const AmgMatrix(5, 5) K = trkCov * R;
245  const AmgMatrix(5, 5) M = s_unitMatrix - K;
246  // --- compute local filtered state
247  const AmgVector(5) newPar = trkPar + K * r;
248  // --- compute filtered covariance matrix
249  const AmgSymMatrix(5) newCov =
250  trkCov.similarity(M) + measCov.similarity(K);
251  const double chiSquared = Amg::chi2(R, r);
252  // create the FQSonSurface
253  fQ = Trk::FitQualityOnSurface(chiSquared, 5);
254  TP.updateParameters(newPar, newCov);
255  return true;
256 }
257 
258 template<int DIM>
259 bool
260 calculateFilterStep_T(Trk::TrackParameters& TP,
261  const AmgSymMatrix(5) & trkCov,
262  const AmgVector(DIM) & measPar,
263  const AmgSymMatrix(DIM) & measCov,
264  int paramKey,
265  Trk::FitQualityOnSurface& fQ)
266 {
267 
268  // get the parameter vector
269  const AmgVector(5)& trkPar = TP.parameters();
270  // reduction matrix
271  const AmgMatrix(DIM, 5) H =
272  s_reMatrices.expansionMatrix(paramKey).topLeftCorner<DIM, 5>();
273  // the projected parameters from the TrackParameters
274  AmgVector(DIM) projTrkPar;
275  if (paramKey == 3 || paramKey == 7 || paramKey == 15) {
276  projTrkPar = trkPar.head<DIM>();
277  } else {
278  projTrkPar = H * trkPar;
279  }
280 
281  // reduction matrix H, Kalman gain K, residual r, combined covariance R
282  // residual after reduction
283  const AmgVector(DIM) r = measPar - projTrkPar;
284  // combined covariance after reduction
285  const AmgSymMatrix(DIM) R =
286  (measCov + projection_T<DIM>(trkCov, paramKey)).inverse();
287  // Kalman gain matrix
288  const AmgMatrix(5, DIM) K = trkCov * H.transpose() * R;
289  const AmgMatrix(5, 5) M = s_unitMatrix - K * H;
290  // --- compute local filtered state
291  const AmgVector(5) newPar = trkPar + K * r;
292  // --- compute filtered covariance matrix
293  // C = M * trkCov * M.T() + K * covRio * K.T()
294  const AmgSymMatrix(5) newCov =
295  trkCov.similarity(M) + K * measCov * K.transpose();
296  const double chiSquared = Amg::chi2(R, r);
297  // create the FQSonSurface
298  fQ = Trk::FitQualityOnSurface(chiSquared, DIM);
299  // In place update of parameters
300  TP.updateParameters(newPar, newCov);
301  return true;
302 }
303 
304 /*
305  * Method implementing the actual kalman filter step
306  * using the methods above.
307  * Also does heavy use of Eigen
308  * so we use flatten
309  */
311 bool
312 filterStep(Trk::TrackParameters& trackParameters,
313  Trk::FitQualityOnSurface& fitQos,
314  const Trk::LocalParameters& measurement,
315  const Amg::MatrixX& measCovariance)
316 {
317 
318  const AmgSymMatrix(5)* trkCov = trackParameters.covariance();
319  if (!trkCov) {
320  return false;
321  }
322 
323  const int nLocCoord = measCovariance.cols();
324  if (!(measurement.dimension() == nLocCoord)) {
325  return false;
326  }
327 
328  switch (nLocCoord) {
329  case 1: {
330  return calculateFilterStep_1D(trackParameters,
331  *trkCov,
332  measurement(0),
333  measCovariance(0, 0),
334  measurement.parameterKey(),
335  fitQos);
336  }
337  case 2: {
338  return calculateFilterStep_T<2>(trackParameters,
339  *trkCov,
340  measurement.head<2>(),
341  measCovariance.topLeftCorner<2, 2>(),
342  measurement.parameterKey(),
343  fitQos);
344  }
345  case 3: {
346  return calculateFilterStep_T<3>(trackParameters,
347  *trkCov,
348  measurement.head<3>(0),
349  measCovariance.topLeftCorner<3, 3>(),
350  measurement.parameterKey(),
351  fitQos);
352  }
353  case 4: {
354  return calculateFilterStep_T<4>(trackParameters,
355  *trkCov,
356  measurement.head<4>(0),
357  measCovariance.topLeftCorner<4, 4>(),
358  measurement.parameterKey(),
359  fitQos);
360  }
361  case 5: {
362  return calculateFilterStep_5D(trackParameters,
363  *trkCov,
364  measurement.head<5>(),
365  measCovariance.topLeftCorner<5, 5>(),
366  fitQos);
367  }
368  default: {
369  return false;
370  }
371  }
372 }
373 
374 /*
375  * Methods needed to calculate the
376  * fit quality
377  */
378 bool
379 makeChi2_1D(Trk::FitQualityOnSurface& updatedFitQoS,
380  const AmgVector(5) & trkPar,
381  const AmgSymMatrix(5) & trkCov,
382  double valRio,
383  double rioCov,
384  int paramKey)
385 {
386 
387  const int mk = measurementCoord_1D(paramKey);
388  const double r = valRio - trkPar(mk);
389  // if (mk==3) catchPiPi;
390  double chiSquared = rioCov + trkCov(mk, mk);
391  if (chiSquared == 0.0) {
392  return false;
393  }
394  chiSquared = r * r / chiSquared;
395  updatedFitQoS = Trk::FitQualityOnSurface(chiSquared, 1);
396  return true;
397 }
398 
399 template<int DIM>
400 bool
401 makeChi2_T(Trk::FitQualityOnSurface& updatedFitQoS,
402  const AmgVector(5) & trkPar,
403  const AmgSymMatrix(5) & trkCov,
404  const AmgVector(DIM) & measPar,
405  const AmgSymMatrix(DIM) & covPar,
406  int paramKey)
407 
408 {
409 
410  const AmgMatrix(DIM, 5) H =
411  s_reMatrices.expansionMatrix(paramKey).topLeftCorner<DIM, 5>();
412  const AmgVector(DIM) r = measPar - H * trkPar;
413  // get the projected matrix
414  AmgSymMatrix(DIM) R = projection_T<DIM>(trkCov, paramKey);
415  R += covPar;
416  // calcualte the chi2 value
417  const double chiSquared = Amg::chi2(R.inverse(), r);
418  updatedFitQoS = Trk::FitQualityOnSurface(chiSquared, DIM);
419  return true;
420 }
421 
422 bool
423 stateFitQuality(Trk::FitQualityOnSurface& updatedFitQoS,
424  const Trk::TrackParameters& trkPar,
425  const Trk::LocalParameters& position,
426  const Amg::MatrixX& covariance)
427 {
428  if (!trkPar.covariance()) {
429  return false;
430  }
431  // For the LocalPos. version, need to get # meas. coord. from covariance
432  // matrix.
433  const int nLocCoord = covariance.cols();
434  switch (nLocCoord) {
435  case 1: {
436  return makeChi2_1D(updatedFitQoS,
437  trkPar.parameters(),
438  (*trkPar.covariance()),
439  position[Trk::locX],
440  covariance(0, 0),
441  1);
442  }
443  case 2: {
444  return makeChi2_T<2>(updatedFitQoS,
445  trkPar.parameters(),
446  (*trkPar.covariance()),
447  position,
448  covariance.topLeftCorner<2, 2>(),
449  3);
450  }
451  default: {
452  return false;
453  }
454  }
455 }
456 
457 /*
458  * Methods to adjust the weights of the
459  * multicomponent state (Gaussian sum)
460  * given a measurement.
461  */
462 template<int DIM>
463 std::pair<double, double>
464 calculateWeight_T(const Trk::TrackParameters* componentTrackParameters,
465  const AmgSymMatrix(5) * predictedCov,
466  const AmgVector(DIM) & measPar,
467  const AmgSymMatrix(DIM) & measCov,
468  int paramKey)
469 {
470  // Define the expansion matrix
471  const AmgMatrix(DIM, 5) H =
472  s_reMatrices.expansionMatrix(paramKey).topLeftCorner<DIM, 5>();
473 
474  // Calculate the residual
475  const AmgVector(DIM) r = measPar - H * componentTrackParameters->parameters();
476  // Residual covariance. Posterior weights is calculated used predicted state
477  // and measurement. Therefore add covariances
478  AmgSymMatrix(DIM) const R(measCov + predictedCov->similarity(H));
479  // compute determinant of residual
480  const double det = R.determinant();
481  if (det == 0) {
482  return { 0, 0 };
483  }
484  // Compute Chi2
485  return {det, (1. / (double)DIM) * Amg::chi2(R.inverse(), r)};
486 }
487 
488 std::pair<double, double>
489 calculateWeight_1D(const Trk::TrackParameters* componentTrackParameters,
490  const AmgSymMatrix(5) * predictedCov,
491  const double measPar,
492  const double measCov,
493  int paramKey)
494 {
495  const int mk = measurementCoord_1D(paramKey);
496  // Calculate the residual
497  const double r = measPar - (componentTrackParameters->parameters())(mk);
498  // Residual covariance. Posterior weights is calculated used predicted state
499  // and measurement. Therefore add covariances
500  const double R = measCov + (*predictedCov)(mk, mk);
501  if (R == 0) {
502  return { 0, 0 };
503  }
504  // Compute Chi2
505  return { R, r * r / R };
506 }
507 
508 std::pair<double, double>
509 calculateWeight_2D_3(const Trk::TrackParameters* componentTrackParameters,
510  const AmgSymMatrix(5) * predictedCov,
511  const AmgVector(2) & measPar,
512  const AmgSymMatrix(2) & measCov)
513 {
514  // Calculate the residual
515  const AmgVector(2) r =
516  measPar - componentTrackParameters->parameters().head<2>();
517  // Residual covariance. Posterior weights is calculated used predicted state
518  // and measurement. Therefore add covariances
519  const AmgSymMatrix(2) R(measCov + predictedCov->topLeftCorner<2, 2>());
520  // compute determinant of residual
521  const double det = R.determinant();
522  if (det == 0) {
523  return { 0, 0 };
524  }
525  // Compute Chi2
526  return {det, 0.5 * Amg::chi2(R.inverse(), r)};
527 }
528 
530 weights(Trk::MultiComponentState&& predictedState,
531  const Trk::MeasurementBase& measurement)
532 {
533 
534  const size_t predictedStateSize = predictedState.size();
535  if (predictedStateSize == 0) {
536  return {};
537  }
538  if (predictedStateSize > GSFConstants::maxNumberofStateComponents) {
539  throw std::runtime_error("PosteriorWeightsCalculator :Invalid predictedState size");
540  }
541  const Trk::LocalParameters& measurementLocalParameters =
542  measurement.localParameters();
543  const int nLocCoord = measurement.localParameters().dimension();
544  if (nLocCoord < 1 || nLocCoord > 5) {
545  return {};
546  }
547 
548  // Move to output and update
549  Trk::MultiComponentState returnMultiComponentState = std::move(predictedState);
550 
551  // Calculate chi2 and determinant of each component.
552  componentsCache determinantRandChi2{};
553  double minimumChi2(10.e10); // Initalise high
554  // Loop over all components
555  for (const auto& component : returnMultiComponentState) {
556 
557  const Trk::TrackParameters* componentTrackParameters =
558  component.params.get();
559  if (!componentTrackParameters) {
560  continue;
561  }
562  const AmgSymMatrix(5)* predictedCov =
563  componentTrackParameters->covariance();
564  if (!predictedCov) {
565  continue;
566  }
567 
568  std::pair<double, double> result(0, 0);
569  switch (nLocCoord) {
570  case 1: {
571  result = calculateWeight_1D(componentTrackParameters,
572  predictedCov,
573  measurementLocalParameters(0),
574  measurement.localCovariance()(0, 0),
575  measurementLocalParameters.parameterKey());
576  } break;
577  case 2: {
578  if (measurementLocalParameters.parameterKey() == 3) {
579  result = calculateWeight_2D_3(
580  componentTrackParameters,
581  predictedCov,
582  measurementLocalParameters.head<2>(),
583  measurement.localCovariance().topLeftCorner<2, 2>());
584  } else {
585  result = calculateWeight_T<2>(
586  componentTrackParameters,
587  predictedCov,
588  measurementLocalParameters.head<2>(),
589  measurement.localCovariance().topLeftCorner<2, 2>(),
590  measurementLocalParameters.parameterKey());
591  }
592  } break;
593  case 3: {
594  result =
595  calculateWeight_T<3>(componentTrackParameters,
596  predictedCov,
597  measurementLocalParameters.head<3>(),
598  measurement.localCovariance().topLeftCorner<3, 3>(),
599  measurementLocalParameters.parameterKey());
600  } break;
601  case 4: {
602  result =
603  calculateWeight_T<4>(componentTrackParameters,
604  predictedCov,
605  measurementLocalParameters.head<4>(),
606  measurement.localCovariance().topLeftCorner<4, 4>(),
607  measurementLocalParameters.parameterKey());
608  } break;
609  case 5: {
610  result =
611  calculateWeight_T<5>(componentTrackParameters,
612  predictedCov,
613  measurementLocalParameters.head<5>(),
614  measurement.localCovariance().topLeftCorner<5, 5>(),
615  measurementLocalParameters.parameterKey());
616  } break;
617  default: {
618  }
619  }
620 
621  if (result.first == 0) {
622  continue;
623  }
624  // Cache R and Chi2
625  determinantRandChi2.elements[determinantRandChi2.numElements] = {
626  result.first, result.second
627  };
628  ++determinantRandChi2.numElements;
629  if (result.second < minimumChi2) {
630  minimumChi2 = result.second;
631  }
632  } // end loop over components
633 
634  // If something went wrong in the loop return empty
635  if (determinantRandChi2.numElements != predictedStateSize) {
636  return {};
637  }
638 
639  // Calculate posterior weights.
640  size_t index(0);
641  double sumWeights(0.);
642  std::array<double, GSFConstants::maxNumberofStateComponents> fallBackWeights{};
643  auto componentItr = returnMultiComponentState.begin();
644  for (; componentItr != returnMultiComponentState.end();
645  ++componentItr, ++index) {
646  // Extract common factor to avoid numerical problems during exponentiation
647  const double chi2 = determinantRandChi2.elements[index].chi2 - minimumChi2;
648  const double priorWeight = componentItr->weight;
649  fallBackWeights[index] = priorWeight;
650  double updatedWeight(0.);
651  // Determinant can not be below 1e-19. Rather ugly but protect
652  // against 0 determinants Normally occur when the component is a poor fit
653  if (determinantRandChi2.elements[index].determinantR > 1e-20) {
654  updatedWeight =
655  priorWeight *
656  sqrt(1. / determinantRandChi2.elements[index].determinantR) *
657  exp(-0.5 * chi2);
658  } else {
659  updatedWeight = 1e-10;
660  }
661  componentItr->weight = updatedWeight;
662  sumWeights += updatedWeight;
663  }
664  if (sumWeights > 0.) {
665  const double invertSumWeights = 1. / sumWeights;
666  // Renormalise the state to total weight = 1
667  for (auto& returnComponent : returnMultiComponentState) {
668  returnComponent.weight *= invertSumWeights;
669  }
670  } else {
671  // If the sum weights is less than 0 revert them back
672  size_t fallbackIndex(0);
673  for (auto& returnComponent : returnMultiComponentState) {
674  returnComponent.weight = fallBackWeights[fallbackIndex];
675  ++fallbackIndex;
676  }
677  }
678  return returnMultiComponentState;
679 }
680 
681 /*
682  * Methods to handle cases with invalid
683  * components
684  */
685 bool
686 invalidComponent(const Trk::TrackParameters* trackParameters)
687 {
688  const auto* measuredCov = trackParameters->covariance();
689  if (!measuredCov) {
690  return true;
691  }
692 
693  for (int i = 0; i < 5; ++i) {
694  if ((*measuredCov)(i, i) <= 0.) {
695  return true;
696  }
697  }
698 
699  return false;
700 }
701 
703 rebuildState(Trk::MultiComponentState&& stateBeforeUpdate)
704 {
705  Trk::MultiComponentState stateWithInsertedErrors =
706  std::move(stateBeforeUpdate);
707  // We need to loop checking for invalid componets i.e negative covariance
708  // diagonal elements and update them with a large covariance matrix
709  for (const Trk::ComponentParameters& component : stateWithInsertedErrors) {
710  const Trk::TrackParameters* trackParameters = component.params.get();
711  const bool rebuildCov = invalidComponent(trackParameters);
712  if (rebuildCov) {
713  AmgSymMatrix(5) bigNewCovarianceMatrix = AmgSymMatrix(5)::Zero();
714  const double covarianceScaler = 1.;
715  bigNewCovarianceMatrix(0, 0) = 250. * covarianceScaler;
716  bigNewCovarianceMatrix(1, 1) = 250. * covarianceScaler;
717  bigNewCovarianceMatrix(2, 2) = 0.25;
718  bigNewCovarianceMatrix(3, 3) = 0.25;
719  bigNewCovarianceMatrix(4, 4) = 0.001 * 0.001;
720  component.params->updateParameters(trackParameters->parameters(),
721  bigNewCovarianceMatrix);
722  }
723  }
724  return stateWithInsertedErrors;
725 }
726 
727 /*
728  * Methods that bring all
729  * weights adjustement, filter step,
730  * and FitQuality together,
731  * so as to update the current
732  * multi-component state given a measurement
733  */
735 calculateFilterStep(Trk::MultiComponentState&& stateBeforeUpdate,
736  const Trk::MeasurementBase& measurement,
737  Trk::FitQualityOnSurface& fitQoS)
738 {
739  // state Assembler cache
741  if (stateBeforeUpdate.empty()) {
742  return {};
743  }
744 
745  // Calculate the weight of each component after the measurement
746  Trk::MultiComponentState stateWithNewWeights =
747  weights(std::move(stateBeforeUpdate), measurement);
748 
749  if (stateWithNewWeights.empty()) {
750  return {};
751  }
752 
753  double chiSquared = 0;
754  int degreesOfFreedom = 0;
755  for (Trk::ComponentParameters& component : stateWithNewWeights) {
756  if (stateWithNewWeights.size() > 1 &&
757  std::abs(component.params->parameters()[Trk::qOverP]) > 0.033333) {
758  continue;
759  }
760  Trk::FitQualityOnSurface componentFitQuality;
762  bool const updateSuccess = filterStep(*(component.params),
763  componentFitQuality,
764  measurement.localParameters(),
765  measurement.localCovariance());
766  if (!updateSuccess) {
767  continue;
768  }
769 
770  if (invalidComponent(component.params.get())) {
771  continue;
772  }
773 
774  if (componentFitQuality.chiSquared() <= 0.) {
775  continue;
776  }
777 
778  chiSquared += component.weight * componentFitQuality.chiSquared();
779 
780  // The same measurement is included in each update
781  // so we can update the degree of freedom only
782  if (degreesOfFreedom == 0.0) {
783  degreesOfFreedom = componentFitQuality.numberDoF();
784  }
785 
786  // Add component to state being prepared for assembly
788  std::move(component));
789  }
790 
791  Trk::MultiComponentState assembledUpdatedState =
793 
794  if (assembledUpdatedState.empty()) {
795  return {};
796  }
797 
798  fitQoS.setChiSquared(chiSquared);
799  fitQoS.setNumberDoF(degreesOfFreedom);
800  // Renormalise state
802  return assembledUpdatedState;
803 }
804 
805 } // end of anonymous namespace
806 
807 // The external/interface methods
808 // implemented in terms of the internal ones.
811  const Trk::MeasurementBase& measurement,
812  FitQualityOnSurface& fitQoS)
813 {
814 
815  // Check all components have associated error matrices
816  bool rebuildStateWithErrors = false;
817  // Perform initial check of state awaiting update. If all states have
818  // associated error matricies then no need to perform the rebuild
819  for (const Trk::ComponentParameters& component : stateBeforeUpdate) {
820  rebuildStateWithErrors =
821  rebuildStateWithErrors || invalidComponent(component.params.get());
822  }
823 
824  if (rebuildStateWithErrors) {
825  Trk::MultiComponentState stateWithInsertedErrors =
826  rebuildState(std::move(stateBeforeUpdate));
827  // Perform the measurement update with the modified state
828  Trk::MultiComponentState updatedState = calculateFilterStep(
829  std::move(stateWithInsertedErrors), measurement, fitQoS);
830  if (updatedState.empty()) {
831  return {};
832  }
833  return updatedState;
834  }
835 
836  // Perform the measurement update
837  Trk::MultiComponentState updatedState = calculateFilterStep(
838  std::move(stateBeforeUpdate), measurement, fitQoS);
839 
840  if (updatedState.empty()) {
841  return {};
842  }
843  return updatedState;
844 }
845 
848  const MeasurementBase& measurement)
849 {
850  double chi2 = 0;
851  Trk::FitQualityOnSurface componentFitQuality;
852  for (const Trk::ComponentParameters& component: updatedState) {
853  stateFitQuality(componentFitQuality,
854  *component.params.get(),
855  measurement.localParameters(),
856  measurement.localCovariance());
857 
858  chi2 += component.weight * componentFitQuality.chiSquared();
859  }
860 
861  //The same measurement is included
862  return { chi2, componentFitQuality.numberDoF() };
863 }
864 
AllowedVariables::e
e
Definition: AsgElectronSelectorTool.cxx:37
beamspotman.r
def r
Definition: beamspotman.py:672
Trk::LocalParameters
Definition: LocalParameters.h:98
Amg::MatrixX
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > MatrixX
Dynamic Matrix - dynamic allocation.
Definition: EventPrimitives.h:27
inline_hints.h
get_generator_info.result
result
Definition: get_generator_info.py:21
MeasurementBase.h
Trk::MultiComponentStateAssembler::assembledState
MultiComponentState assembledState(MultiComponentStateAssembler::Cache &&cache)
Method to return the cached state object - it performs a reweighting before returning the object base...
Definition: MultiComponentStateAssembler.cxx:120
Trk::differentialCheck
@ differentialCheck
Definition: KalmanUpdatorSMatrix.h:58
Trk::locX
@ locX
Definition: ParamDefs.h:37
Trk::LocalParameters::parameterKey
int parameterKey() const
Identifier key for matrix expansion/reduction.
TRTCalib_Extractor.det
det
Definition: TRTCalib_Extractor.py:36
index
Definition: index.py:1
Trk::ProjectionMatricesSet
the matrices to access the variably-dimensioned local parameters and map them to the defined five tra...
Definition: ProjectionMatricesSet.h:29
theta
Scalar theta() const
theta method
Definition: AmgMatrixBasePlugin.h:75
taskman.template
dictionary template
Definition: taskman.py:314
M_PI
#define M_PI
Definition: ActiveFraction.h:11
Trk::MultiComponentStateAssembler::Cache
Definition: MultiComponentStateAssembler.h:35
const
bool const RAWDATA *ch2 const
Definition: LArRodBlockPhysicsV0.cxx:560
drawFromPickle.exp
exp
Definition: drawFromPickle.py:36
Trk::ParametersCommon::updateParameters
void updateParameters(const AmgVector(DIM) &, const AmgSymMatrix(DIM) &)
Update parameters and covariance , passing covariance by ref.
AmgSymMatrix
#define AmgSymMatrix(dim)
Definition: EventPrimitives.h:50
AmgMatrix
#define AmgMatrix(rows, cols)
Definition: EventPrimitives.h:49
similarity
Matrix< Scalar, OtherDerived::RowsAtCompileTime, OtherDerived::RowsAtCompileTime > similarity(const MatrixBase< OtherDerived > &m) const
similarity method : yields ms = m*s*m^T
Definition: AmgMatrixBasePlugin.h:133
xAOD::TrackParameters
TrackParameters_v1 TrackParameters
Definition: Event/xAOD/xAODTracking/xAODTracking/TrackParameters.h:11
GSFConstants::maxNumberofStateComponents
constexpr int8_t maxNumberofStateComponents
The state is described by N Gaussian components The Beth Heitler Material effect are also described b...
Definition: GsfConstants.h:43
Trk::FitQualityOnSurface
Definition: FitQualityOnSurface.h:19
MultiComponentStateAssembler.h
H
#define H(x, y, z)
Definition: MD5.cxx:114
lumiFormat.i
int i
Definition: lumiFormat.py:85
ATH_FLATTEN
#define ATH_FLATTEN
Definition: inline_hints.h:52
Trk::theta
@ theta
Definition: ParamDefs.h:66
AmgVector
AmgVector(4) T2BSTrackFilterTool
Definition: T2BSTrackFilterTool.cxx:114
chi2
double chi2(TH1 *h0, TH1 *h1)
Definition: comparitor.cxx:525
Trk::MultiComponentState
std::vector< ComponentParameters > MultiComponentState
Definition: ComponentParameters.h:27
AnalysisUtils::Delta::R
double R(const INavigable4Momentum *p1, const double v_eta, const double v_phi)
Definition: AnalysisMisc.h:49
Trk::ParametersBase
Definition: ParametersBase.h:55
xAOD::double
double
Definition: CompositeParticle_v1.cxx:159
Trk::FitQualityOnSurface::numberDoF
int numberDoF() const
returns the number of degrees of freedom of the overall track or vertex fit as integer
Definition: FitQuality.h:60
Trk::RangeCheckDef
RangeCheckDef
Definition: KalmanUpdatorSMatrix.h:56
Trk::MeasurementBase::localCovariance
const Amg::MatrixX & localCovariance() const
Interface method to get the localError.
Definition: MeasurementBase.h:138
Trk::MeasurementBase
Definition: MeasurementBase.h:58
head
std::string head(std::string s, const std::string &pattern)
head of a string
Definition: computils.cxx:312
Trk
Ensure that the ATLAS eigen extensions are properly loaded.
Definition: FakeTrackBuilder.h:9
Amg
Definition of ATLAS Math & Geometry primitives (Amg)
Definition: AmgStringHelpers.h:19
Trk::MultiComponentStateHelpers::renormaliseState
void renormaliseState(MultiComponentState &, double norm=1)
Performing renormalisation of total state weighting to one.
Definition: ComponentParameters.cxx:80
TRT::Track::degreesOfFreedom
@ degreesOfFreedom
Definition: InnerDetector/InDetCalibEvent/TRT_CalibData/TRT_CalibData/TrackInfo.h:79
Trk::absoluteCheck
@ absoluteCheck
Definition: KalmanUpdatorSMatrix.h:57
weights
Definition: herwig7_interface.h:38
Trk::ComponentParameters
Definition: ComponentParameters.h:22
LocalParameters.h
Trk::MeasurementBase::localParameters
const LocalParameters & localParameters() const
Interface method to get the LocalParameters.
Definition: MeasurementBase.h:132
IDTPM::chiSquared
float chiSquared(const U &p)
Definition: TrackParametersHelper.h:128
Trk::GsfMeasurementUpdator::fitQuality
FitQualityOnSurface fitQuality(const MultiComponentState &, const MeasurementBase &)
Method for determining the chi2 of the multi-component state and the number of degrees of freedom.
Definition: GsfMeasurementUpdator.cxx:847
Trk::LocalParameters::dimension
int dimension() const
Dimension of this localParameters() vector.
EventPrimitivesCovarianceHelpers.h
DeMoScan.index
string index
Definition: DeMoScan.py:362
Trk::qOverP
@ qOverP
perigee
Definition: ParamDefs.h:67
if
if(febId1==febId2)
Definition: LArRodBlockPhysicsV0.cxx:567
physics_parameters.parameters
parameters
Definition: physics_parameters.py:144
GsfConstants.h
Trk::phi
@ phi
Definition: ParamDefs.h:75
Trk::GsfMeasurementUpdator::update
MultiComponentState update(Trk::MultiComponentState &&, const Trk::MeasurementBase &, FitQualityOnSurface &fitQoS)
Method for updating the multi-state with a new measurement and calculate the fit qaulity at the same ...
Definition: GsfMeasurementUpdator.cxx:810
PowhegControl_ttFCNC_NLO.params
params
Definition: PowhegControl_ttFCNC_NLO.py:226
FitQuality.h
GsfMeasurementUpdator.h
Code for performing kalman filter update step on multi-component states for the gaussian-sum filter.
Trk::MultiComponentStateAssembler::addComponent
void addComponent(MultiComponentStateAssembler::Cache &cache, ComponentParameters &&multiComponentState)
Method to add a single set of Trk::ComponentParameters to the cached Trk::MultiComponentState object ...
Definition: MultiComponentStateAssembler.h:71
Trk::FitQualityOnSurface::chiSquared
double chiSquared() const
returns the of the overall track fit
Definition: FitQuality.h:56
fitman.k
k
Definition: fitman.py:528
generate::Zero
void Zero(TH1D *hin)
Definition: generate.cxx:32
Amg::chi2
double chi2(const T &precision, const U &residual, const int sign=1)
Definition: EventPrimitivesCovarianceHelpers.h:221
mapkey::key
key
Definition: TElectronEfficiencyCorrectionTool.cxx:37