Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2026-04-06 16:37:22

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/EffPlotTool.hpp"
0010 
0011 #include "Acts/Utilities/VectorHelpers.hpp"
0012 #include "Acts/Utilities/Zip.hpp"
0013 #include "ActsExamples/EventData/SimParticle.hpp"
0014 
0015 #include <format>
0016 #include <limits>
0017 
0018 using Acts::VectorHelpers::eta;
0019 using Acts::VectorHelpers::perp;
0020 using Acts::VectorHelpers::phi;
0021 
0022 namespace ActsExamples {
0023 
0024 EffPlotTool::EffPlotTool(const EffPlotTool::Config& cfg,
0025                          Acts::Logging::Level lvl)
0026     : m_cfg(cfg), m_logger(Acts::getDefaultLogger("EffPlotTool", lvl)) {
0027   ACTS_DEBUG("Initialize the histograms for efficiency plots");
0028 
0029   // 1D efficiencies
0030   m_efficiencies1D.insert(
0031       {"trackeff_vs_eta",
0032        Efficiency1("trackeff_vs_eta",
0033                    std::format("Tracking efficiency with pT > {} GeV/c",
0034                                m_cfg.minTruthPt / Acts::UnitConstants::GeV),
0035                    std::array{m_cfg.varBinning.at("Eta")})});
0036   m_efficiencies1D.insert(
0037       {"trackeff_vs_phi",
0038        Efficiency1("trackeff_vs_phi",
0039                    std::format("Tracking efficiency with pT > {} GeV/c",
0040                                m_cfg.minTruthPt / Acts::UnitConstants::GeV),
0041                    std::array{m_cfg.varBinning.at("Phi")})});
0042   m_efficiencies1D.insert(
0043       {"trackeff_vs_pT", Efficiency1("trackeff_vs_pT", "Tracking efficiency",
0044                                      std::array{m_cfg.varBinning.at("Pt")})});
0045   m_efficiencies1D.insert(
0046       {"trackeff_vs_LogPt",
0047        Efficiency1("trackeff_vs_LogPt", "Tracking efficiency",
0048                    std::array{m_cfg.varBinning.at("LogPt")})});
0049   m_efficiencies1D.insert(
0050       {"trackeff_vs_LowPt",
0051        Efficiency1("trackeff_vs_LowPt", "Tracking efficiency",
0052                    std::array{m_cfg.varBinning.at("LowPt")})});
0053   m_efficiencies1D.insert(
0054       {"trackeff_vs_d0",
0055        Efficiency1("trackeff_vs_d0",
0056                    std::format("Tracking efficiency with pT > {} GeV/c",
0057                                m_cfg.minTruthPt / Acts::UnitConstants::GeV),
0058                    std::array{m_cfg.varBinning.at("D0")})});
0059   m_efficiencies1D.insert(
0060       {"trackeff_vs_z0",
0061        Efficiency1("trackeff_vs_z0",
0062                    std::format("Tracking efficiency with pT > {} GeV/c",
0063                                m_cfg.minTruthPt / Acts::UnitConstants::GeV),
0064                    std::array{m_cfg.varBinning.at("Z0")})});
0065   m_efficiencies1D.insert(
0066       {"trackeff_vs_DeltaR",
0067        Efficiency1("trackeff_vs_DeltaR",
0068                    std::format("Tracking efficiency with pT > {} GeV/c",
0069                                m_cfg.minTruthPt / Acts::UnitConstants::GeV),
0070                    std::array{m_cfg.varBinning.at("DeltaR")})});
0071   m_efficiencies1D.insert(
0072       {"trackeff_vs_prodR",
0073        Efficiency1("trackeff_vs_prodR",
0074                    std::format("Tracking efficiency with pT > {} GeV/c",
0075                                m_cfg.minTruthPt / Acts::UnitConstants::GeV),
0076                    std::array{m_cfg.varBinning.at("prodR")})});
0077 
0078   // 2D efficiencies
0079   m_efficiencies2D.insert(
0080       {"trackeff_vs_eta_phi",
0081        Efficiency2("trackeff_vs_eta_phi",
0082                    std::format("Tracking efficiency with pT > {} GeV/c",
0083                                m_cfg.minTruthPt / Acts::UnitConstants::GeV),
0084                    std::array{m_cfg.varBinning.at("Eta"),
0085                               m_cfg.varBinning.at("Phi")})});
0086   m_efficiencies2D.insert(
0087       {"trackeff_vs_eta_pt",
0088        Efficiency2(
0089            "trackeff_vs_eta_pt", "Tracking efficiency",
0090            std::array{m_cfg.varBinning.at("Eta"), m_cfg.varBinning.at("Pt")})});
0091 
0092   const auto& etaAxis = m_cfg.varBinning.at("Eta");
0093   const auto& ptAxis = m_cfg.varBinning.at("Pt");
0094 
0095   // efficiency vs eta in different pT ranges
0096   for (const auto& [i, ptRange] : Acts::enumerate(m_cfg.truthPtRangesForEta)) {
0097     const std::string name = std::format("trackeff_vs_eta_ptRange_{}", i);
0098     const std::string title =
0099         std::format("Tracking efficiency with pT in [{}, {}] GeV/c",
0100                     ptRange.first / Acts::UnitConstants::GeV,
0101                     ptRange.second / Acts::UnitConstants::GeV);
0102     m_trackEffVsEtaInPtRanges.emplace_back(name, title, std::array{etaAxis});
0103   }
0104 
0105   // efficiency vs pT in different abs(eta) ranges
0106   for (const auto& [i, absEtaRange] :
0107        Acts::enumerate(m_cfg.truthAbsEtaRangesForPt)) {
0108     const std::string name = std::format("trackeff_vs_pT_absEtaRange_{}", i);
0109     const std::string title =
0110         std::format("Tracking efficiency with |#eta| in [{}, {}]",
0111                     absEtaRange.first, absEtaRange.second);
0112     m_trackEffVsPtInAbsEtaRanges.emplace_back(name, title, std::array{ptAxis});
0113   }
0114 }
0115 
0116 void EffPlotTool::fill(const Acts::GeometryContext& gctx,
0117                        const SimParticleState& truthParticle,
0118                        const double deltaR, const bool status) {
0119   constexpr double nan = std::numeric_limits<double>::quiet_NaN();
0120 
0121   const auto intersection =
0122       m_cfg.beamline
0123           ->intersect(gctx, truthParticle.position(), truthParticle.direction())
0124           .closest();
0125   Acts::Vector2 d0z0{nan, nan};
0126   if (intersection.isValid()) {
0127     auto localRes = m_cfg.beamline->globalToLocal(gctx, intersection.position(),
0128                                                   truthParticle.direction());
0129     if (localRes.ok()) {
0130       d0z0 = localRes.value();
0131     }
0132   }
0133 
0134   const double t_phi = phi(truthParticle.direction());
0135   const double t_eta = eta(truthParticle.direction());
0136   const double t_absEta = std::abs(t_eta);
0137   const double t_pT = truthParticle.transverseMomentum();
0138   const double t_d0 = d0z0.x();
0139   const double t_z0 = d0z0.y();
0140   const double t_deltaR = deltaR;
0141   const double t_prodR = perp(truthParticle.position());
0142 
0143   // cut on truth pT with the global range for the relevant plots
0144   if (t_pT >= m_cfg.minTruthPt) {
0145     m_efficiencies1D.at("trackeff_vs_eta").fill({t_eta}, status);
0146     m_efficiencies1D.at("trackeff_vs_phi").fill({t_phi}, status);
0147     m_efficiencies1D.at("trackeff_vs_d0").fill({t_d0}, status);
0148     m_efficiencies1D.at("trackeff_vs_z0").fill({t_z0}, status);
0149     m_efficiencies1D.at("trackeff_vs_DeltaR").fill({t_deltaR}, status);
0150     m_efficiencies1D.at("trackeff_vs_prodR").fill({t_prodR}, status);
0151 
0152     m_efficiencies2D.at("trackeff_vs_eta_phi").fill({t_eta, t_phi}, status);
0153   }
0154 
0155   // do not cut on truth pT as it is a variable on the plot
0156   m_efficiencies1D.at("trackeff_vs_pT").fill({t_pT}, status);
0157   m_efficiencies1D.at("trackeff_vs_LogPt").fill({t_pT}, status);
0158   m_efficiencies1D.at("trackeff_vs_LowPt").fill({t_pT}, status);
0159   m_efficiencies2D.at("trackeff_vs_eta_pt").fill({t_eta, t_pT}, status);
0160 
0161   // fill the efficiency vs eta in different pT ranges
0162   for (auto&& [ptRange, eff] :
0163        Acts::zip(m_cfg.truthPtRangesForEta, m_trackEffVsEtaInPtRanges)) {
0164     if (t_pT >= ptRange.first && t_pT < ptRange.second) {
0165       eff.fill({t_eta}, status);
0166     }
0167   }
0168 
0169   // fill the efficiency vs pT in different eta ranges
0170   for (auto&& [absEtaRange, eff] :
0171        Acts::zip(m_cfg.truthAbsEtaRangesForPt, m_trackEffVsPtInAbsEtaRanges)) {
0172     if (t_absEta >= absEtaRange.first && t_absEta < absEtaRange.second) {
0173       eff.fill({t_pT}, status);
0174     }
0175   }
0176 }
0177 
0178 }  // namespace ActsExamples