ATLAS Offline Software
GaussianSumFitter.cxx
Go to the documentation of this file.
1 /*ยง
2  Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3 */
4 
13 //
16 //
27 #include "TrkTrack/Track.h"
28 #include "TrkTrack/TrackInfo.h"
29 //
30 #include <algorithm>
31 #include <vector>
32 
33 namespace {
41 std::unique_ptr<Trk::FitQuality> buildFitQuality(
42  const Trk::GaussianSumFitter::GSFTrajectory& smoothedTrajectory) {
43  double chiSquared = 0.;
44  int numberDoF = -5;
45  Trk::GaussianSumFitter::GSFTrajectory::const_iterator stateOnSurface =
46  smoothedTrajectory.begin();
47  for (; stateOnSurface != smoothedTrajectory.end(); ++stateOnSurface) {
48  if (!stateOnSurface->typeFlags.test(Trk::TrackStateOnSurface::Measurement)) {
49  continue;
50  }
51  if (!stateOnSurface->fitQualityOnSurface) {
52  continue;
53  }
54  chiSquared += stateOnSurface->fitQualityOnSurface.chiSquared();
55  numberDoF += stateOnSurface->fitQualityOnSurface.numberDoF();
56  }
57  if (std::isnan(chiSquared) || chiSquared <= 0.) {
58  return nullptr;
59  }
60  return std::make_unique<Trk::FitQuality>(chiSquared, numberDoF);
61 }
62 
67 GSFTsos smootherHelper(Trk::MultiComponentState&& updatedState,
68  std::unique_ptr<Trk::MeasurementBase>&& measurement,
70  bool combineToSingle, bool useMode) {
71  if (combineToSingle) {
73  updatedState, useMode);
74  if (combinedLastState) {
75  return {fitQuality, std::move(measurement), std::move(combinedLastState),
76  std::move(updatedState)};
77  }
78  }
79  return {fitQuality, std::move(measurement), nullptr, std::move(updatedState)};
80 }
81 
87 Trk::MeasurementSet stripMeasurements(
88  const Trk::Track& inputTrack, const Trk::MeasurementSet* addMeasurementSet,
89  bool reintegrateOutliers) {
90  Trk::MeasurementSet measurementSet;
91  for (const auto* tsos : *(inputTrack.trackStateOnSurfaces())) {
92  if (!tsos) {
93  continue;
94  }
95  bool use =
97  (reintegrateOutliers && tsos->type(Trk::TrackStateOnSurface::Outlier));
98  if (!use) {
99  continue;
100  }
101  const Trk::MeasurementBase* meas = tsos->measurementOnTrack();
102  if (meas) {
103  measurementSet.push_back(meas);
104  }
105  }
106 
107  if (addMeasurementSet) {
108  for (const auto* meas : *addMeasurementSet) {
109  measurementSet.push_back(meas);
110  }
111  }
112  return measurementSet;
113 }
114 
120 Trk::PrepRawDataSet stripPrepRawData(
121  const Trk::Track& inputTrack, const Trk::PrepRawDataSet* addPrepRawDataSet,
122  bool reintegrateOutliers) {
123 
124  Trk::PrepRawDataSet prepRawDataSet;
125  for (const auto* tsos : *(inputTrack.trackStateOnSurfaces())) {
126  if (!tsos) {
127  continue;
128  }
129  bool use =
131  (reintegrateOutliers && tsos->type(Trk::TrackStateOnSurface::Outlier));
132  if (!use) {
133  continue;
134  }
135  const Trk::MeasurementBase* meas = tsos->measurementOnTrack();
136  if (!meas) {
137  continue;
138  }
139  const Trk::RIO_OnTrack* rioOnTrack = nullptr;
141  rioOnTrack = static_cast<const Trk::RIO_OnTrack*>(meas);
142  }
143  if (!rioOnTrack) {
144  continue;
145  }
146  const Trk::PrepRawData* prepRawData = rioOnTrack->prepRawData();
147  if (prepRawData) {
148  prepRawDataSet.push_back(prepRawData);
149  }
150  }
151  if (addPrepRawDataSet) {
152  for (const auto* prepRawData : *addPrepRawDataSet) {
153  prepRawDataSet.push_back(prepRawData);
154  }
155  }
156  return prepRawDataSet;
157 }
158 
159 } // end of anonymous namespace
160 
162  const std::string& name,
163  const IInterface* parent)
165  , m_trkParametersComparisonFunction{}
166 {
167  declareInterface<ITrackFitter>(this);
168 }
169 
172 {
173 
174  if (m_maximumNumberOfComponents > GSFConstants::maxNumberofStateComponents) {
175  ATH_MSG_FATAL("Requested MaximumNumberOfComponents > "
177  return StatusCode::FAILURE;
178  }
179  // GSF extrapolator
180  ATH_CHECK(m_extrapolator.retrieve());
181  // We need a to create RIO_OnTrack (calibrated/corrected)
182  // measurements only if we start from PrePRawData.
183  if (!m_refitOnMeasurementBase) {
184  ATH_CHECK(m_rioOnTrackCreator.retrieve());
185  } else {
186  m_rioOnTrackCreator.disable();
187  }
188  // Initialise the closest track parameters search algorithm
189  const Amg::Vector3D referencePosition(0, 0, 0);
190  m_trkParametersComparisonFunction = Trk::TrkParametersComparisonFunction(referencePosition);
191  m_particleHypothesis = m_extrapolator->particleHypothesis();
192  return StatusCode::SUCCESS;
193 }
194 
199 std::unique_ptr<Trk::Track>
201  const EventContext& ctx,
202  const Trk::Track& inputTrack,
203  const Trk::RunOutlierRemoval outlierRemoval,
204  const Trk::ParticleHypothesis particleHypothesis) const
205 {
206  // Check that the input track has well defined parameters
207  if (inputTrack.trackParameters()->empty()) {
208  ATH_MSG_FATAL("No track parameters on Track!");
209  return nullptr;
210  }
211  // Check that the input track has associated MeasurementBase objects
212  if (inputTrack.trackStateOnSurfaces()->empty()) {
213  ATH_MSG_FATAL("Empty MeasurementBase ");
214  return nullptr;
215  }
216  // Retrieve the set of track parameters closest to the reference point
217  const Trk::TrackParameters* paramNearestRef =
218  *(std::min_element(inputTrack.trackParameters()->begin(),
219  inputTrack.trackParameters()->end(),
220  m_trkParametersComparisonFunction));
221 
222  // Rrefit using measurement base then
223  if (m_refitOnMeasurementBase) {
224  MeasurementSet measurementSet =
225  stripMeasurements(inputTrack, nullptr, m_reintegrateOutliers);
226  return fit(ctx, measurementSet, *paramNearestRef, outlierRemoval,
227  particleHypothesis);
228  }
229  // Refit using PrepRawData level then
230  PrepRawDataSet prepRawDataSet =
231  stripPrepRawData(inputTrack, nullptr, m_reintegrateOutliers);
232  return fit(ctx, prepRawDataSet, *paramNearestRef, outlierRemoval,
233  particleHypothesis);
234 }
235 
239 std::unique_ptr<Trk::Track>
240 Trk::GaussianSumFitter::fit(const EventContext& ctx,
241  const Track& intrk,
242  const PrepRawDataSet& addPrdColl,
243  const RunOutlierRemoval runOutlier,
244  const ParticleHypothesis matEffects) const
245 {
246  // Check that the input track has well defined parameters
247  if (intrk.trackParameters()->empty()) {
248  ATH_MSG_FATAL("No track parameters on Track!");
249  return nullptr;
250  }
251  // Check that the input track has associated MeasurementBase objects
252  if (intrk.trackStateOnSurfaces()->empty()) {
253  ATH_MSG_FATAL("Empty MeasurementBase ");
254  return nullptr;
255  }
256  // determine the Track Parameter closest to the reference
257  const TrackParameters* paramNearestRef = *(std::min_element(
258  intrk.trackParameters()->begin(), intrk.trackParameters()->end(),
259  m_trkParametersComparisonFunction));
260  PrepRawDataSet prepRawDataSet =
261  stripPrepRawData(intrk, &addPrdColl, m_reintegrateOutliers);
262  return fit(ctx, prepRawDataSet, *paramNearestRef, runOutlier, matEffects);
263 }
264 
268 std::unique_ptr<Trk::Track>
269 Trk::GaussianSumFitter::fit(const EventContext& ctx,
270  const Track& inputTrack,
271  const MeasurementSet& addSet,
272  const RunOutlierRemoval runOutlier,
273  const ParticleHypothesis matEffects) const
274 {
275  // Check that the input track has well defined parameters
276  if (inputTrack.trackParameters()->empty()) {
277  ATH_MSG_FATAL("No estimation of track parameters near origin!");
278  return nullptr;
279  }
280  // Check that the input track has associated MeasurementBase objects
281  if (inputTrack.trackStateOnSurfaces()->empty()) {
283  "Attempting to fit track to empty MeasurementBase "
284  "collection!");
285  return nullptr;
286  }
287  // Retrieve the set of track parameters closest to the reference point
288  const Trk::TrackParameters* paramNearestRef = *(std::min_element(
289  inputTrack.trackParameters()->begin(),
290  inputTrack.trackParameters()->end(), m_trkParametersComparisonFunction));
291 
292  MeasurementSet measurementSet =
293  stripMeasurements(inputTrack, &addSet, m_reintegrateOutliers);
294  return fit(ctx, measurementSet, *paramNearestRef, runOutlier, matEffects);
295 }
296 
301 std::unique_ptr<Trk::Track> Trk::GaussianSumFitter::fit(
302  const EventContext& ctx,
303  const Track& intrk1,
304  const Track& intrk2,
305  const RunOutlierRemoval runOutlier,
306  const ParticleHypothesis matEffects) const {
307  // protection against not having measurements on the input tracks
308  if (!intrk1.trackStateOnSurfaces() || !intrk2.trackStateOnSurfaces() ||
309  intrk1.trackStateOnSurfaces()->size() < 2) {
311  "called to refit empty track or track with too little "
312  "information, reject fit");
313  return nullptr;
314  }
315  if (!intrk1.trackParameters() || intrk1.trackParameters()->empty()) {
317  "input #1 fails to provide track parameters for "
318  "seeding the GXF, reject fit");
319  return nullptr;
320  }
321  // Here we assume that the reference paramrs are the beginning
322  // of the 1st track.
323  const TrackParameters* minPar = *(intrk1.trackParameters()->begin());
324  MeasurementSet measurementSet1 =
325  stripMeasurements(intrk1, nullptr, m_reintegrateOutliers);
326  MeasurementSet measurementSet2 =
327  stripMeasurements(intrk2, nullptr, m_reintegrateOutliers);
328 
329  for (const auto* meas : measurementSet2) {
330  measurementSet1.push_back(meas);
331  }
332 
333  return fit(ctx, measurementSet1, *minPar, runOutlier, matEffects);
334 }
335 
340 std::unique_ptr<Trk::Track>
342  const EventContext& ctx,
343  const Trk::PrepRawDataSet& prepRawDataSet,
344  const Trk::TrackParameters& estimatedParametersNearOrigin,
345  const Trk::RunOutlierRemoval /* Not used*/,
346  const Trk::ParticleHypothesis /* Not used*/) const
347 {
348  if (prepRawDataSet.size() < 3) {
349  ATH_MSG_WARNING("Requesting Fit with less than three prep Raw Data!!!");
350  return nullptr;
351  }
352  // We need a sorted PrepRawDataSet
353  Trk::PrepRawDataSet sortedPrepRawDataSet = PrepRawDataSet(prepRawDataSet);
354  const Trk::PrepRawDataComparisonFunction prdComparisonFunction =
356  estimatedParametersNearOrigin.position(),
357  estimatedParametersNearOrigin.momentum());
358 
359  std::sort(sortedPrepRawDataSet.begin(), sortedPrepRawDataSet.end(),
360  prdComparisonFunction);
361 
362  // Create Extrapolator cache that holds material effects cache;
363  Trk::IMultiStateExtrapolator::Cache extrapolatorCache;
364  // Perform GSF forward fit
365  GSFTrajectory forwardTrajectory =
366  forwardFit(ctx, extrapolatorCache, sortedPrepRawDataSet,
367  estimatedParametersNearOrigin);
368 
369  if (forwardTrajectory.empty()) {
370  return nullptr;
371  }
372  // Perform GSF smoother operation
373  GSFTrajectory smoothedTrajectory = smootherFit(
374  ctx, extrapolatorCache, forwardTrajectory);
375  if (smoothedTrajectory.empty()) {
376  return nullptr;
377  }
378 
379  // Fit quality
380  std::unique_ptr<FitQuality> fitQuality = buildFitQuality(smoothedTrajectory);
381  if (!fitQuality) {
382  return nullptr;
383  }
384 
385  // Create parameters at perigee if needed
386  auto perigeeMultiStateOnSurface = makePerigee(
387  ctx, extrapolatorCache, smoothedTrajectory);
388  if (!perigeeMultiStateOnSurface.multiComponentState.empty()) {
389  smoothedTrajectory.push_back(std::move(perigeeMultiStateOnSurface));
390  } else {
391  return nullptr;
392  }
393 
394  // Reverse the order of the TSOS's after the smoother
395  // to make be order flow from inside to out
396  std::reverse(smoothedTrajectory.begin(), smoothedTrajectory.end());
397 
398  // Create Trk::Track
400  info.setTrackProperties(TrackInfo::BremFit);
401  info.setTrackProperties(TrackInfo::BremFitSuccessful);
402  return std::make_unique<Track>(info, convertTrajToTrack(smoothedTrajectory),
403  std::move(fitQuality));
404 }
405 
410 std::unique_ptr<Trk::Track>
412  const EventContext& ctx,
413  const Trk::MeasurementSet& measurementSet,
414  const Trk::TrackParameters& estimatedParametersNearOrigin,
415  const Trk::RunOutlierRemoval /* Not used*/,
416  const Trk::ParticleHypothesis /*Not used*/) const
417 {
418  if (measurementSet.size() < 3) {
419  ATH_MSG_WARNING("Requesting fit with less than 3 Measurements!!!");
420  return nullptr;
421  }
422  // We need to separate the possible CaloCluster on Track
423  // measurement from the rest
424  const Trk::CaloCluster_OnTrack* ccot(nullptr);
425 
426  Trk::MeasurementSet cleanedMeasurementSet;
427  cleanedMeasurementSet.reserve(measurementSet.size());
428 
429  for (const auto* meas : measurementSet) {
430  if (!meas) {
431  ATH_MSG_WARNING("There is an empty MeasurementBase object ! ");
432  continue;
433  }
435  ccot = static_cast<const Trk::CaloCluster_OnTrack*>(meas);
436  } else {
437  cleanedMeasurementSet.push_back(meas);
438  }
439  }
440 
441  // We need a sorted measurement set
442  Trk::MeasurementSet sortedMeasurementSet =
443  MeasurementSet(cleanedMeasurementSet);
444 
446  measurementBaseComparisonFunction(
447  estimatedParametersNearOrigin.position(),
448  estimatedParametersNearOrigin.momentum());
449  sort(sortedMeasurementSet.begin(), sortedMeasurementSet.end(),
450  measurementBaseComparisonFunction);
451 
452  // Create Extrapolator cache that holds material effects cache;
453  Trk::IMultiStateExtrapolator::Cache extrapolatorCache;
454 
455  // Perform GSF forwards fit
456  GSFTrajectory forwardTrajectory =
457  forwardFit(ctx, extrapolatorCache, sortedMeasurementSet,
458  estimatedParametersNearOrigin);
459 
460  if (forwardTrajectory.empty()) {
461  return nullptr;
462  }
463 
464  // Perform GSF smoother operation
465  GSFTrajectory smoothedTrajectory = smootherFit(
466  ctx, extrapolatorCache, forwardTrajectory, ccot);
467  if (smoothedTrajectory.empty()) {
468  return nullptr;
469  }
470 
471  // fit quality
472  std::unique_ptr<FitQuality> fitQuality = buildFitQuality(smoothedTrajectory);
473  if (!fitQuality) {
474  return nullptr;
475  }
476 
477  // Create parameters at perigee if needed
478  auto perigeeMultiStateOnSurface = makePerigee(
479  ctx, extrapolatorCache, smoothedTrajectory);
480  if (!perigeeMultiStateOnSurface.multiComponentState.empty()) {
481  smoothedTrajectory.push_back(std::move(perigeeMultiStateOnSurface));
482  } else {
483  return nullptr;
484  }
485 
486  // Reverse the order of the TSOS's to make be order flow from inside to out
487  std::reverse(smoothedTrajectory.begin(), smoothedTrajectory.end());
488 
489  // Create track
491  info.setTrackProperties(TrackInfo::BremFit);
492  info.setTrackProperties(TrackInfo::BremFitSuccessful);
493  return std::make_unique<Track>(info, convertTrajToTrack(smoothedTrajectory),
494  std::move(fitQuality));
495 }
496 
500 std::unique_ptr<MultiComponentStateOnSurfaceDV> Trk::GaussianSumFitter::convertTrajToTrack(
501  GSFTrajectory& trajectory) const{
502  const bool slimTransientMTSOS = m_slimTransientMTSOS;
503  auto MTSOS = std::make_unique<MultiComponentStateOnSurfaceDV>();
504  MTSOS->reserve(trajectory.size());
505  for (GSFTsos& state : trajectory) {
506  MTSOS->push_back(state.convert(slimTransientMTSOS));
507  }
508  return MTSOS;
509 }
510 
514 GSFTsos
516  const EventContext& ctx,
517  Trk::IMultiStateExtrapolator::Cache& extrapolatorCache,
518  const GSFTrajectory& smoothedTrajectory) const
519 {
520 
521  // Start at the end of the smoothed trajectory
522  // we should be the closest to perigee/origin.
523  const GSFTsos&
524  multiComponentStateOnSurfaceNearestOrigin = smoothedTrajectory.back();
525  const Trk::MultiComponentState* multiComponentState =
526  &(multiComponentStateOnSurfaceNearestOrigin.multiComponentState);
527 
528  // Extrapolate to perigee
529  const Trk::PerigeeSurface perigeeSurface;
530  Trk::MultiComponentState stateExtrapolatedToPerigee =
531  m_extrapolator->extrapolate(ctx,
532  extrapolatorCache,
533  *multiComponentState,
534  perigeeSurface,
536  false);
537 
538  if (stateExtrapolatedToPerigee.empty()) {
539  return {};
540  }
541 
542  // Determine the combined state
543  std::unique_ptr<Trk::TrackParameters> combinedPerigee =
544  MultiComponentStateCombiner::combineToSingle(stateExtrapolatedToPerigee, m_useMode);
545 
546  // Perigee is an additional MultiComponentStateOnSurface
547  std::bitset<Trk::TrackStateOnSurface::NumberOfTrackStateOnSurfaceTypes>
548  pattern(0);
550 
551  if(std::abs(combinedPerigee->position().z())>5000.) {
552  ATH_MSG_WARNING("Pathological perigee well outside of tracking detector!! Returning {}");
553  return {};
554  }
555 
556  if (std::abs(combinedPerigee->parameters()[Trk::qOverP]) > 1e8) {
558  "makePerigee() about to return with 0 momentum!! Returning {}");
559  return {};
560  }
561 
562  return {Trk::FitQualityOnSurface{},
563  nullptr,
564  std::move(combinedPerigee),
565  std::move(stateExtrapolatedToPerigee),
566  pattern};
567 }
568 
569 /*
570  * Implementation of the smoothing of the trajectory.
571  * 1. We start from the forward trajectory containing
572  * the measurement and prediction from the fwd fit.
573  * 2. We do an update for the last state of the forward fit,
574  * as this is special for the EDM
575  * 3. Then we inflate the last state covariance matrix and we
576  * run a filter in the backward direction.
577  * 4. We optionally can use a two filter "Bayessian" smoother
578  * ,see "Track fitting with non-Gaussian noise" R.Fruhwirth,
579  * since we kept the "predictedStates" in the forward pass.
580  *
581  * This is also our chance to do any special treatement
582  * needed for the ATLAS EDM.
583  */
586  const EventContext& ctx,
587  Trk::IMultiStateExtrapolator::Cache& extrapolatorCache,
588  GSFTrajectory& forwardTrajectory,
589  const Trk::CaloCluster_OnTrack* ccot) const
590 {
591  if (forwardTrajectory.empty()) {
593  "Attempting to smooth an empty forward trajectory... Exiting!");
594  return {};
595  }
596 
597  GSFTrajectory smoothedTrajectory;
598  smoothedTrajectory.reserve(forwardTrajectory.size());
599  // For the smoother we start from the end and we go
600  // back. We need to find the the first track
601  // state on surface in this reverse direction
602  auto trackStateOnSurfaceItr = forwardTrajectory.rbegin();
603  bool foundMeasurement = false;
604  for (; trackStateOnSurfaceItr != forwardTrajectory.rend(); ++trackStateOnSurfaceItr) {
605  if (!(*trackStateOnSurfaceItr).typeFlags.test(TrackStateOnSurface::Measurement)) {
606  smoothedTrajectory.emplace_back(std::move(*trackStateOnSurfaceItr));
607  } else {
608  foundMeasurement = true;
609  break;
610  }
611  }
612  if(!foundMeasurement){
613  return {};
614  }
615  // This is the 1st track state on surface for a measurement
616  // in the reverse direction. Our starting point for the smoother
617  auto smootherPredictionMultiState = std::move(trackStateOnSurfaceItr->multiComponentState);
618  // Perform an update with the measurement.
619  // This will be the first entry in the trajectory
620  std::unique_ptr<Trk::MeasurementBase> firstSmootherMeasurementOnTrack =
621  std::move(trackStateOnSurfaceItr->measurementOnTrack);
622  if (!firstSmootherMeasurementOnTrack) {
624  "Initial state on surface in smoother does not have an associated "
625  "MeasurementBase object");
626  return {};
627  }
629  Trk::MultiComponentState firstSmoothedState =
630  Trk::GsfMeasurementUpdator::update(std::move(smootherPredictionMultiState),
631  *firstSmootherMeasurementOnTrack,
632  fitQuality);
633  if (firstSmoothedState.empty()) {
634  return {};
635  }
636  if (!MultiComponentStateHelpers::allHaveCovariance(firstSmoothedState)) {
638  "Not all components have covariance. Rejecting smoothed state.");
639  return {};
640  }
641  // The first in reverse (last in normal order) TSOS is also special
642  // for the EDM. Sso we also a proper collapse
643  // of the multi component to single TrackParameter
644  std::unique_ptr<Trk::TrackParameters> combinedFirstSmoothedState =
646  m_useMode);
647  smoothedTrajectory.emplace_back(
648  fitQuality, std::move(firstSmootherMeasurementOnTrack),
649  std::move(combinedFirstSmoothedState),
650  MultiComponentStateHelpers::clone(firstSmoothedState));
651  const auto& updatedFirstStateOnSurface = smoothedTrajectory.back();
652  //
653  // Now we are ready to proceed
654  // Generate a prediction with an inflated covariance of all components
655  // in the first smoothed state.
656  // This way there is no dependance on error of prediction
657  //
658  // NB local Y and theta need care due to TRT.
659  //
660  Trk::MultiComponentState smoothedStateWithScaledError =
662  std::move(firstSmoothedState), m_smootherCovFactors[0],
663  m_smootherCovFactors[1], m_smootherCovFactors[2],
664  m_smootherCovFactors[3], m_smootherCovFactors[4]);
665  // Do the 1st update using the state with inflated covariance.
666  Trk::FitQualityOnSurface fitQualityWithScaledErrors;
668  std::move(smoothedStateWithScaledError),
669  *(updatedFirstStateOnSurface.measurementOnTrack),
670  fitQualityWithScaledErrors);
671  if (updatedState.empty()) {
672  ATH_MSG_WARNING("Smoother prediction could not be determined");
673  return {};
674  }
675  // loopUpdatedState points to the most recent predicted MultiComponentState.
676  Trk::MultiComponentState* loopUpdatedState = &updatedState;
677  // Increase reverse iterator by one.
678  ++trackStateOnSurfaceItr;
679  // The is the last one we will see as we go back.
680  auto lasttrackStateOnSurface = forwardTrajectory.rend() - 1;
681  // And this is the pre-last
682  auto secondLastTrackStateOnSurface = forwardTrajectory.rend() - 2;
683  //Loop
684  for (; trackStateOnSurfaceItr != forwardTrajectory.rend();++trackStateOnSurfaceItr) {
685  auto& trackStateOnSurface = (*trackStateOnSurfaceItr);
686  // Retrieve the MeasurementBase object from the TrackStateOnSurface object
687  auto measurement = std::move(trackStateOnSurface.measurementOnTrack);
688  if (!measurement) {
689  ATH_MSG_WARNING("MeasurementBase object could not be extracted from a "
690  "measurement TSOS...continuing");
691  continue;
692  }
693  // Create prediction for the next measurement surface.
694  // For the smoother the direction of propagation
695  // is opposite to the direction of momentum
696  Trk::MultiComponentState extrapolatedState = m_extrapolator->extrapolate(
697  ctx, extrapolatorCache, (*loopUpdatedState),
698  measurement->associatedSurface(), Trk::oppositeMomentum, false);
699 
700  if (extrapolatedState.empty()) {
701  return {};
702  }
703  // Handle the case where Original measurement was flagged as an outlier.
704  if (!trackStateOnSurface.typeFlags.test(TrackStateOnSurface::Measurement)) {
705  std::bitset<TrackStateOnSurface::NumberOfTrackStateOnSurfaceTypes> type(
706  0);
708  smoothedTrajectory.emplace_back(FitQualityOnSurface(1, 1),
709  std::move(measurement), nullptr,
710  std::move(extrapolatedState), type);
711  loopUpdatedState = &(smoothedTrajectory.back().multiComponentState);
712  continue;
713  }
714  // Update with the measurement
715  updatedState = Trk::GsfMeasurementUpdator::update(
716  std::move(extrapolatedState), *measurement, fitQuality);
717  if (updatedState.empty()) {
718  ATH_MSG_WARNING("Could not update the multi-component state");
719  return {};
720  }
721  // last in reverse (first in normal order) is special for the EDM
722  const bool islast = (trackStateOnSurfaceItr == lasttrackStateOnSurface);
723  if (m_combineWithFitter) {
724  // Optional combine smoother state with fitter state
725  const Trk::MultiComponentState& forwardsMultiState =
726  trackStateOnSurface.multiComponentState;
727  Trk::MultiComponentState combinedfitterState =
728  Trk::TwoStateCombiner::combine(forwardsMultiState, updatedState,
729  m_maximumNumberOfComponents);
730  if (combinedfitterState.empty()) {
732  "Could not combine state from forward fit with "
733  "smoother state");
734  return {};
735  }
736  auto combinedFitQuality = Trk::GsfMeasurementUpdator::fitQuality(
737  combinedfitterState, *measurement);
738  smoothedTrajectory.emplace_back(
739  smootherHelper(std::move(combinedfitterState), std::move(measurement),
740  combinedFitQuality, islast, m_useMode));
741  } else {
742  // If combination with forwards state is not done
743  smoothedTrajectory.emplace_back(
744  smootherHelper(std::move(updatedState), std::move(measurement),
745  fitQuality, islast, m_useMode));
746  }
747  // Handle adding measurement from calo if it is present
748  if (ccot && trackStateOnSurfaceItr == secondLastTrackStateOnSurface) {
749  if (!addCCOT(ctx, ccot, smoothedTrajectory)) {
750  ATH_MSG_WARNING("Could not add Calo Cluster On Track Measurement");
751  return {};
752  }
753  }
754  // For the next iteration start from last added
755  loopUpdatedState = &(smoothedTrajectory.back().multiComponentState);
756  } // End for loop over all components
757  return smoothedTrajectory;
758 }
759 
764 bool
766  const EventContext& ctx,
767  const Trk::CaloCluster_OnTrack* ccot,
768  GSFTrajectory& smoothedTrajectory) const
769 {
770  const GSFTsos& currentMultiStateOS = smoothedTrajectory.back();
771  if (!ccot) {
772  return false;
773  }
774  const auto& currentMultiComponentState = currentMultiStateOS.multiComponentState;
775  const Trk::MeasurementBase* measurement = currentMultiStateOS.measurementOnTrack.get();
776  if (!measurement) {
777  return false;
778  }
779  const Trk::Surface& currentSurface = measurement->associatedSurface();
780 
781  //Create a ccot that we will own. So we use it to build our own TSOS
782  auto ownCCOT = ccot->uniqueClone();
783  // Extrapolate to the Calo to get prediction
784  Trk::MultiComponentState extrapolatedToCaloState = m_extrapolator->extrapolateDirectly(
785  ctx, currentMultiComponentState, ownCCOT->associatedSurface(),
786  Trk::alongMomentum, false);
787 
788  if (extrapolatedToCaloState.empty()) {
789  return false;
790  }
791  // Update newly extrapolated state with measurement
793  Trk::MultiComponentState updatedStateAtCalo =
794  Trk::GsfMeasurementUpdator::update(std::move(extrapolatedToCaloState),
795  *ownCCOT, fitQuality);
796  if (updatedStateAtCalo.empty()) {
797  return false;
798  }
799 
800  // Extrapolate back to the surface near the origin
801  auto improvedState = m_extrapolator->extrapolateDirectly(
802  ctx, updatedStateAtCalo, currentSurface, Trk::oppositeMomentum,
803  false);
804 
805  if (improvedState.empty()) {
806  return false;
807  }
808  // Combine the improved state after extrapolating back from the calo
809  // and find the mode of the distribution
810  std::unique_ptr<Trk::TrackParameters> combinedSingleState =
811  MultiComponentStateCombiner::combineToSingle(improvedState, m_useMode);
812 
813  // Now build a dummy measurement for the improved estimation
815  covMatrix.setZero();
816  covMatrix(0, 0) = 1e6;
818  Trk::LocalParameters locpars(locX);
819  auto pseudoMeasurement = std::make_unique<Trk::PseudoMeasurementOnTrack>(
820  std::move(locpars), std::move(covMatrix), currentSurface);
821 
822  // Build TSOS at the surface of calo
823  smoothedTrajectory.emplace_back(
824  fitQuality,
825  std::move(ownCCOT),
826  nullptr,
827  std::move(updatedStateAtCalo));
828 
829  // Build a TSOS using the dummy measurement and and the final combined state
830  smoothedTrajectory.emplace_back(
832  // We do not add the fitquality again. As this would be the one we
833  // add at calo. Here not a real measurement
834  std::move(pseudoMeasurement),
835  std::move(combinedSingleState),
836  std::move(improvedState));
837 
838  return true;
839 }
mergePhysValFiles.pattern
pattern
Definition: DataQuality/DataQualityUtils/scripts/mergePhysValFiles.py:25
Trk::MultiComponentStateHelpers::allHaveCovariance
bool allHaveCovariance(const MultiComponentState &in)
Check to see if all components have covariance Matrix.
Definition: ComponentParameters.cxx:66
Trk::LocalParameters
Definition: LocalParameters.h:98
EstimatedBremOnTrack.h
TwoStateCombiner.h
Trk::TrackInfo
Contains information about the 'fitter' of this track.
Definition: Tracking/TrkEvent/TrkTrack/TrkTrack/TrackInfo.h:32
Trk::TrackStateOnSurface::Perigee
@ Perigee
This represents a perigee, and so will contain a Perigee object only.
Definition: TrackStateOnSurface.h:117
ATH_MSG_FATAL
#define ATH_MSG_FATAL(x)
Definition: AthMsgStreamMacros.h:34
Trk::TrackInfo::GaussianSumFilter
@ GaussianSumFilter
Tracks from Gaussian Sum Filter.
Definition: Tracking/TrkEvent/TrkTrack/TrkTrack/TrackInfo.h:53
Trk::PrepRawDataSet
std::vector< const PrepRawData * > PrepRawDataSet
vector of clusters and drift circles
Definition: FitterTypes.h:26
TrackParameters.h
PerigeeSurface.h
Trk::locX
@ locX
Definition: ParamDefs.h:37
GaussianSumFitter.h
Class for fitting according to the Gaussian Sum Filter Formalism Implements Algorithm C described in ...
Trk::Track
The ATLAS Track class.
Definition: Tracking/TrkEvent/TrkTrack/TrkTrack/Track.h:73
Trk::GaussianSumFitter::initialize
virtual StatusCode initialize() override final
AlgTool initialise method.
Definition: GaussianSumFitter.cxx:171
Trk::PerigeeSurface
Definition: PerigeeSurface.h:43
Trk::ParametersBase::position
const Amg::Vector3D & position() const
Access method for the position.
Trk::oppositeMomentum
@ oppositeMomentum
Definition: PropDirection.h:21
Trk::Track::trackStateOnSurfaces
const Trk::TrackStates * trackStateOnSurfaces() const
return a pointer to a const DataVector of const TrackStateOnSurfaces.
Trk::MeasurementBaseType::CaloCluster_OnTrack
@ CaloCluster_OnTrack
Definition: MeasurementBase.h:53
Trk::TrackInfo::BremFitSuccessful
@ BremFitSuccessful
A brem fit was performed on this track and this fit was successful.
Definition: Tracking/TrkEvent/TrkTrack/TrkTrack/TrackInfo.h:81
Trk::GaussianSumFitter::convertTrajToTrack
std::unique_ptr< MultiComponentStateOnSurfaceDV > convertTrajToTrack(GSFTrajectory &trajectory) const
Helper to convert the GSFTrajectory to a Trk::Track.
Definition: GaussianSumFitter.cxx:500
Trk::TwoStateCombiner::combine
Trk::MultiComponentState combine(const Trk::MultiComponentState &forwardsMultiState, const Trk::MultiComponentState &smootherMultiState, unsigned int maximumNumberOfComponents)
Helper to combine forward with smoother MultiComponentStates.
Definition: TwoStateCombiner.cxx:149
Trk::RIO_OnTrack
Definition: RIO_OnTrack.h:70
Trk::GaussianSumFitter::smootherFit
GSFTrajectory smootherFit(const EventContext &ctx, Trk::IMultiStateExtrapolator::Cache &, GSFTrajectory &forwardTrajectory, const CaloCluster_OnTrack *ccot=nullptr) const
Gsf smoothed trajectory.
Definition: GaussianSumFitter.cxx:585
Trk::alongMomentum
@ alongMomentum
Definition: PropDirection.h:20
Trk::MultiComponentStateCombiner::combineToSingle
std::unique_ptr< Trk::TrackParameters > combineToSingle(const MultiComponentState &, const bool useMode=false)
@bried Calculate combined state of many components
Definition: MultiComponentStateCombiner.cxx:190
PrepRawData.h
python.CaloAddPedShiftConfig.type
type
Definition: CaloAddPedShiftConfig.py:42
Trk::RunOutlierRemoval
bool RunOutlierRemoval
switch to toggle quality processing after fit
Definition: FitterTypes.h:22
Trk::DefinedParameter
std::pair< double, ParamDefs > DefinedParameter
Definition: DefinedParameter.h:27
Trk::GaussianSumFitter::GSFTrajectory
std::vector< GSFTsos > GSFTrajectory
Definition: GaussianSumFitter.h:124
DeMoUpdate.reverse
reverse
Definition: DeMoUpdate.py:563
Trk::TrackStateOnSurface::Outlier
@ Outlier
This TSoS contains an outlier, that is, it contains a MeasurementBase/RIO_OnTrack which was not used ...
Definition: TrackStateOnSurface.h:122
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
Track.h
Trk::GaussianSumFitter::makePerigee
GSFTsos makePerigee(const EventContext &ctx, Trk::IMultiStateExtrapolator::Cache &, const GSFTrajectory &smoothedTrajectory) const
Produces a perigee from a smoothed trajectory.
Definition: GaussianSumFitter.cxx:515
Trk::GaussianSumFitter::fit
virtual std::unique_ptr< Track > fit(const EventContext &ctx, const PrepRawDataSet &, const TrackParameters &, const RunOutlierRemoval, const ParticleHypothesis) const override final
Fit a collection of 'PrepRawData' objects using the Gaussian Sum Filter.
Definition: GaussianSumFitter.cxx:341
Trk::AmgSymMatrix
AmgSymMatrix(5) &GXFTrackState
Definition: GXFTrackState.h:156
Trk::FitQualityOnSurface
Definition: FitQualityOnSurface.h:19
Trk::ParticleHypothesis
ParticleHypothesis
Definition: ParticleHypothesis.h:28
ATH_MSG_ERROR
#define ATH_MSG_ERROR(x)
Definition: AthMsgStreamMacros.h:33
Trk::CaloCluster_OnTrack::uniqueClone
std::unique_ptr< CaloCluster_OnTrack > uniqueClone() const
NVI Clone giving up unique pointer.
Definition: CaloCluster_OnTrack.h:58
EL::StatusCode
::StatusCode StatusCode
StatusCode definition for legacy code.
Definition: PhysicsAnalysis/D3PDTools/EventLoop/EventLoop/StatusCode.h:22
xAOD::covMatrix
covMatrix
Definition: TrackMeasurement_v1.cxx:19
PrepRawDataComparisonFunction.h
GSFTsos::measurementOnTrack
std::unique_ptr< Trk::MeasurementBase > measurementOnTrack
Definition: GSFTsos.h:24
Trk::MultiComponentStateHelpers::clone
MultiComponentState clone(const MultiComponentState &in)
Clone TrackParameters method.
Definition: ComponentParameters.cxx:15
PseudoMeasurementOnTrack.h
test_pyathena.parent
parent
Definition: test_pyathena.py:15
Trk::GaussianSumFitter::GaussianSumFitter
GaussianSumFitter(const std::string &, const std::string &, const IInterface *)
Constructor with parameters to be passed to AlgTool.
Definition: GaussianSumFitter.cxx:161
Trk::MultiComponentState
std::vector< ComponentParameters > MultiComponentState
Definition: ComponentParameters.h:27
ATH_CHECK
#define ATH_CHECK
Definition: AthCheckMacros.h:40
Trk::MeasurementBase::type
virtual bool type(MeasurementBaseType::Type type) const =0
Interface method checking the type.
Trk::TrackInfo::BremFit
@ BremFit
A brem fit was performed on this track.
Definition: Tracking/TrkEvent/TrkTrack/TrkTrack/TrackInfo.h:78
Trk::ParametersBase
Definition: ParametersBase.h:55
Trk::IMultiStateExtrapolator::Cache
MultiStateExtrapolator cache class.
Definition: IMultiStateExtrapolator.h:51
MuonValidation_CreateResolutionProfiles.fit
def fit(h, emin, emax)
Definition: MuonValidation_CreateResolutionProfiles.py:69
GSFTsos
Definition: GSFTsos.h:17
Trk::MultiComponentStateHelpers::WithScaledError
MultiComponentState WithScaledError(MultiComponentState &&in, double errorScaleLocX, double errorScaleLocY, double errorScalePhi, double errorScaleTheta, double errorScaleQoverP)
Scale the covariance matrix components by individual factors.
Definition: ComponentParameters.cxx:26
Trk::MeasurementSet
std::vector< const MeasurementBase * > MeasurementSet
vector of fittable measurements
Definition: FitterTypes.h:30
Trk::PrepRawData
Definition: PrepRawData.h:62
Trk::MeasurementBase
Definition: MeasurementBase.h:58
Trk::Track::trackParameters
const DataVector< const TrackParameters > * trackParameters() const
Return a pointer to a vector of TrackParameters.
Definition: Tracking/TrkEvent/TrkTrack/src/Track.cxx:97
Trk::CaloCluster_OnTrack
Definition: CaloCluster_OnTrack.h:32
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:240
Trk::MeasurementBase::associatedSurface
virtual const Surface & associatedSurface() const =0
Interface method to get the associated Surface.
TrackInfo.h
RIO_OnTrack.h
Trk::GaussianSumFitter::addCCOT
bool addCCOT(const EventContext &ctx, const Trk::CaloCluster_OnTrack *ccot, GSFTrajectory &smoothedTrajectory) const
Methof to add the CaloCluster onto the track.
Definition: GaussianSumFitter.cxx:765
Trk::pseudoMeasurement
@ pseudoMeasurement
Definition: MeasurementType.h:26
GSFTsos::multiComponentState
Trk::MultiComponentState multiComponentState
Definition: GSFTsos.h:20
Amg::Vector3D
Eigen::Matrix< double, 3, 1 > Vector3D
Definition: GeoPrimitives.h:47
Trk::RIO_OnTrack::prepRawData
virtual const Trk::PrepRawData * prepRawData() const =0
returns the PrepRawData (also known as RIO) object to which this RIO_OnTrack is associated.
Trk::MeasurementBaseType::RIO_OnTrack
@ RIO_OnTrack
Definition: MeasurementBase.h:49
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:842
Trk::ParametersBase::momentum
const Amg::Vector3D & momentum() const
Access method for the momentum.
CaloCluster_OnTrack.h
ATH_MSG_WARNING
#define ATH_MSG_WARNING(x)
Definition: AthMsgStreamMacros.h:32
Trk::TrkParametersComparisonFunction
ComparisonFunction< TrackParameters > TrkParametersComparisonFunction
Definition: TrkParametersComparisonFunction.h:26
Trk::qOverP
@ qOverP
perigee
Definition: ParamDefs.h:67
Track
Definition: TriggerChamberClusterOnTrackCreator.h:21
GsfConstants.h
Trk::PrepRawDataComparisonFunction
Class providing comparison function, or relational definition, for PrepRawData.
Definition: PrepRawDataComparisonFunction.h:33
Trk::MeasurementBaseComparisonFunction
Class implementing a comparison function for sorting MeasurementBase objects.
Definition: MeasurementBaseComparisonFunction.h:40
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:805
AthAlgTool
Definition: AthAlgTool.h:26
FitQuality.h
Trk::Surface
Definition: Tracking/TrkDetDescr/TrkSurfaces/TrkSurfaces/Surface.h:79
DataVector::size
size_type size() const noexcept
Returns the number of elements in the collection.
python.ParticleTypeUtil.info
def info
Definition: ParticleTypeUtil.py:87
DataVector::empty
bool empty() const noexcept
Returns true if the collection is empty.
Trk::TrackStateOnSurface::Measurement
@ Measurement
This is a measurement, and will at least contain a Trk::MeasurementBase.
Definition: TrackStateOnSurface.h:101
MeasurementBaseComparisonFunction.h