Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2026-04-18 07:36:08

0001 // This file is part of the ACTS project.
0002 //
0003 // Copyright (C) 2016 CERN for the benefit of the ACTS project
0004 //
0005 // This Source Code Form is subject to the terms of the Mozilla Public
0006 // License, v. 2.0. If a copy of the MPL was not distributed with this
0007 // file, You can obtain one at https://mozilla.org/MPL/2.0/.
0008 
0009 #include "ActsExamples/Validation/TrackFinderPerformanceCollector.hpp"
0010 
0011 #include "Acts/EventData/BoundTrackParameters.hpp"
0012 #include "Acts/Utilities/Logger.hpp"
0013 #include "Acts/Utilities/VectorHelpers.hpp"
0014 
0015 #include <utility>
0016 
0017 namespace ActsExamples {
0018 
0019 TrackFinderPerformanceCollector::TrackFinderPerformanceCollector(
0020     Config cfg, std::unique_ptr<const Acts::Logger> logger)
0021     : m_cfg(std::move(cfg)),
0022       m_logger(std::move(logger)),
0023       m_effPlotTool(m_cfg.effPlotToolConfig, m_logger->level()),
0024       m_fakePlotTool(m_cfg.fakePlotToolConfig, m_logger->level()),
0025       m_duplicationPlotTool(m_cfg.duplicationPlotToolConfig, m_logger->level()),
0026       m_trackSummaryPlotTool(m_cfg.trackSummaryPlotToolConfig,
0027                              m_logger->level()),
0028       m_trackQualityPlotTool(m_cfg.trackQualityPlotToolConfig,
0029                              m_logger->level()) {
0030   for (const auto& [key, _] : m_cfg.subDetectorTrackSummaryVolumes) {
0031     TrackSummaryPlotTool::Config subConfig = m_cfg.trackSummaryPlotToolConfig;
0032     subConfig.prefix = key;
0033     m_subDetectorSummaryTools.emplace(
0034         std::piecewise_construct, std::forward_as_tuple(key),
0035         std::forward_as_tuple(subConfig, m_logger->level()));
0036   }
0037 }
0038 
0039 ProcessCode TrackFinderPerformanceCollector::fill(
0040     const Acts::GeometryContext& geoContext, const ConstTrackContainer& tracks,
0041     const SimParticleContainer& particles,
0042     const TrackParticleMatching& trackParticleMatching,
0043     const ParticleTrackMatching& particleTrackMatching,
0044     const InverseMultimap<SimBarcode>& particleMeasurementsMap) {
0045   std::size_t unmatched = 0;
0046   std::size_t missingRefSurface = 0;
0047 
0048   for (const auto& track : tracks) {
0049     m_stats.nTotalTracks++;
0050 
0051     if (!track.hasReferenceSurface()) {
0052       missingRefSurface++;
0053       continue;
0054     }
0055 
0056     Acts::BoundTrackParameters fittedParameters =
0057         track.createParametersAtReference();
0058 
0059     m_trackSummaryPlotTool.fill(fittedParameters, track.nTrackStates(),
0060                                 track.nMeasurements(), track.nOutliers(),
0061                                 track.nHoles(), track.nSharedHits());
0062 
0063     for (const auto& [key, volumes] : m_cfg.subDetectorTrackSummaryVolumes) {
0064       std::size_t nTrackStates{};
0065       std::size_t nMeasurements{};
0066       std::size_t nOutliers{};
0067       std::size_t nHoles{};
0068       std::size_t nSharedHits{};
0069 
0070       for (auto state : track.trackStatesReversed()) {
0071         if (!state.hasReferenceSurface() ||
0072             !volumes.contains(state.referenceSurface().geometryId().volume())) {
0073           continue;
0074         }
0075         nTrackStates++;
0076         nMeasurements +=
0077             static_cast<std::size_t>(state.typeFlags().isMeasurement());
0078         nOutliers += static_cast<std::size_t>(state.typeFlags().isOutlier());
0079         nHoles += static_cast<std::size_t>(state.typeFlags().isHole());
0080         nSharedHits +=
0081             static_cast<std::size_t>(state.typeFlags().isSharedHit());
0082       }
0083       m_subDetectorSummaryTools.at(key).fill(fittedParameters, nTrackStates,
0084                                              nMeasurements, nOutliers, nHoles,
0085                                              nSharedHits);
0086     }
0087 
0088     auto imatched = trackParticleMatching.find(track.index());
0089     if (imatched == trackParticleMatching.end()) {
0090       unmatched++;
0091       continue;
0092     }
0093 
0094     const auto& particleMatch = imatched->second;
0095 
0096     if (particleMatch.classification == TrackMatchClassification::Fake) {
0097       m_stats.nTotalFakeTracks++;
0098     }
0099     if (particleMatch.classification == TrackMatchClassification::Duplicate) {
0100       m_stats.nTotalDuplicateTracks++;
0101     }
0102 
0103     m_fakePlotTool.fill(fittedParameters, particleMatch.classification ==
0104                                               TrackMatchClassification::Fake);
0105     m_duplicationPlotTool.fill(
0106         fittedParameters,
0107         particleMatch.classification == TrackMatchClassification::Duplicate);
0108 
0109     if (particleMatch.particle.has_value() &&
0110         particleMeasurementsMap.contains(particleMatch.particle.value())) {
0111       const auto measurements =
0112           particleMeasurementsMap.equal_range(particleMatch.particle.value());
0113 
0114       std::size_t nTrackMeasurements =
0115           track.nMeasurements() + track.nOutliers();
0116       std::size_t nMatchedHits =
0117           particleMatch.contributingParticles.front().hitCount;
0118       std::size_t nParticleHits =
0119           std::distance(measurements.first, measurements.second);
0120 
0121       double completeness = static_cast<double>(nMatchedHits) / nParticleHits;
0122       double purity = static_cast<double>(nMatchedHits) / nTrackMeasurements;
0123 
0124       m_trackQualityPlotTool.fill(fittedParameters, completeness, purity);
0125     }
0126   }
0127 
0128   if (unmatched > 0) {
0129     ACTS_VERBOSE("No matching information found for " << unmatched
0130                                                       << " tracks");
0131   }
0132   if (missingRefSurface > 0) {
0133     ACTS_VERBOSE("Reference surface was missing for " << missingRefSurface
0134                                                       << " tracks");
0135   }
0136 
0137   for (const auto& particle : particles) {
0138     auto particleId = particle.particleId();
0139 
0140     std::size_t nMatchedTracks = 0;
0141     std::size_t nFakeTracks = 0;
0142     bool isReconstructed = false;
0143     if (auto imatched = particleTrackMatching.find(particleId);
0144         imatched != particleTrackMatching.end()) {
0145       isReconstructed = imatched->second.track.has_value();
0146       nMatchedTracks = (isReconstructed ? 1 : 0) + imatched->second.duplicates;
0147 
0148       m_stats.nTotalMatchedTracks += nMatchedTracks;
0149       m_stats.nTotalMatchedParticles += isReconstructed ? 1 : 0;
0150 
0151       if (nMatchedTracks > 1) {
0152         m_stats.nTotalDuplicateParticles += 1;
0153       }
0154 
0155       nFakeTracks = imatched->second.fakes;
0156       if (nFakeTracks > 0) {
0157         m_stats.nTotalFakeParticles += 1;
0158       }
0159     }
0160 
0161     double minDeltaR = -1;
0162     for (const auto& closeParticle : particles) {
0163       if (closeParticle.particleId() == particleId) {
0164         continue;
0165       }
0166       double distance = Acts::VectorHelpers::deltaR(particle.direction(),
0167                                                     closeParticle.direction());
0168       if (minDeltaR == -1 || distance < minDeltaR) {
0169         minDeltaR = distance;
0170       }
0171     }
0172 
0173     m_effPlotTool.fill(geoContext, particle.initialState(), minDeltaR,
0174                        isReconstructed);
0175     m_duplicationPlotTool.fill(particle.initialState(), nMatchedTracks);
0176     m_fakePlotTool.fill(particle.initialState(), nMatchedTracks, nFakeTracks);
0177 
0178     m_stats.nTotalParticles += 1;
0179   }
0180 
0181   return ProcessCode::SUCCESS;
0182 }
0183 
0184 void TrackFinderPerformanceCollector::logSummary() const {
0185   const Acts::Logger& log = *m_logger;
0186   float eff_tracks =
0187       static_cast<float>(m_stats.nTotalMatchedTracks) / m_stats.nTotalTracks;
0188   float fakeRatio_tracks =
0189       static_cast<float>(m_stats.nTotalFakeTracks) / m_stats.nTotalTracks;
0190   float duplicationRatio_tracks =
0191       static_cast<float>(m_stats.nTotalDuplicateTracks) / m_stats.nTotalTracks;
0192 
0193   float eff_particle = static_cast<float>(m_stats.nTotalMatchedParticles) /
0194                        m_stats.nTotalParticles;
0195   float fakeRatio_particle =
0196       static_cast<float>(m_stats.nTotalFakeParticles) / m_stats.nTotalParticles;
0197   float duplicationRatio_particle =
0198       static_cast<float>(m_stats.nTotalDuplicateParticles) /
0199       m_stats.nTotalParticles;
0200 
0201   ACTS_LOG_WITH_LOGGER(
0202       log, Acts::Logging::DEBUG,
0203       "nTotalTracks                = " << m_stats.nTotalTracks);
0204   ACTS_LOG_WITH_LOGGER(
0205       log, Acts::Logging::DEBUG,
0206       "nTotalMatchedTracks         = " << m_stats.nTotalMatchedTracks);
0207   ACTS_LOG_WITH_LOGGER(
0208       log, Acts::Logging::DEBUG,
0209       "nTotalDuplicateTracks       = " << m_stats.nTotalDuplicateTracks);
0210   ACTS_LOG_WITH_LOGGER(
0211       log, Acts::Logging::DEBUG,
0212       "nTotalFakeTracks            = " << m_stats.nTotalFakeTracks);
0213 
0214   ACTS_LOG_WITH_LOGGER(
0215       log, Acts::Logging::INFO,
0216       "Efficiency with tracks (nMatchedTracks/ nAllTracks) = " << eff_tracks);
0217   ACTS_LOG_WITH_LOGGER(
0218       log, Acts::Logging::INFO,
0219       "Fake ratio with tracks (nFakeTracks/nAllTracks) = " << fakeRatio_tracks);
0220   ACTS_LOG_WITH_LOGGER(
0221       log, Acts::Logging::INFO,
0222       "Duplicate ratio with tracks (nDuplicateTracks/nAllTracks) = "
0223           << duplicationRatio_tracks);
0224   ACTS_LOG_WITH_LOGGER(
0225       log, Acts::Logging::INFO,
0226       "Efficiency with particles (nMatchedParticles/nTrueParticles) = "
0227           << eff_particle);
0228   ACTS_LOG_WITH_LOGGER(
0229       log, Acts::Logging::INFO,
0230       "Fake ratio with particles (nFakeParticles/nTrueParticles) = "
0231           << fakeRatio_particle);
0232   ACTS_LOG_WITH_LOGGER(
0233       log, Acts::Logging::INFO,
0234       "Duplicate ratio with particles (nDuplicateParticles/nTrueParticles) = "
0235           << duplicationRatio_particle);
0236 }
0237 
0238 }  // namespace ActsExamples