ATLAS Offline Software
Loading...
Searching...
No Matches
TrackTruthMatchingBaseAlg.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
5#ifndef ACTSTRK_TRACKTRUTHMATCHINGBASEALG_H
6#define ACTSTRK_TRACKTRUTHMATCHINGBASEALG_H 1
7
8// Base Class
10
11// Gaudi includes
12#include "Gaudi/Property.h"
13
14// Handle Keys
17
22
23#include <mutex>
26
27#include <string>
28#include <memory>
29#include <array>
30#include <atomic>
31#include <type_traits>
32
33#include <cmath>
34#include <iomanip>
35#include <ostream>
36#include <sstream>
37#include <vector>
38
39namespace ActsTrk
40{
41 constexpr bool TrackFindingValidationDebugHists = false;
43
45 {
46 template <bool DetailEnabled>
47 struct BaseStat;
48
49 template <bool DetailEnabled>
50 friend struct BaseStat;
51
52 public:
53 TrackTruthMatchingBaseAlg(const std::string &name,
54 ISvcLocator *pSvcLocator);
55
56 virtual StatusCode initialize() override;
57 virtual StatusCode finalize() override;
58
59 protected:
60 const TruthParticleHitCounts &getTruthParticleHitCounts(const EventContext &ctx) const {
61 SG::ReadHandle<TruthParticleHitCounts> truth_particle_hit_counts_handle = SG::makeHandle(m_truthHitCounts, ctx);
62 if (!truth_particle_hit_counts_handle.isValid()) {
63 ATH_MSG_ERROR("No truth particle hit count map for key " << m_truthHitCounts.key() );
64 std::runtime_error("Failed to get truth particle hit count map");
65 }
66 return *truth_particle_hit_counts_handle;
67 }
68
70
71 std::size_t perEtaSize() const {
72 return m_detailedStat.perEtaSize();
73 }
74 std::size_t perPdgIdSize() const {
75 return m_detailedStat.perPdgIdSize();
76 }
77
78 template <bool DetailEnabled>
79 struct EventStatBase : public BaseStat<DetailEnabled> {
80 static constexpr bool doDetail = DetailEnabled;
81
82 EventStatBase(const IAthSelectionTool &truth_selection_tool,
83 std::size_t per_eta_size,
84 std::size_t per_pdg_size,
85 [[maybe_unused]] std::size_t track_to_truth_size)
86 : BaseStat<DetailEnabled>(truth_selection_tool, per_eta_size, per_pdg_size),
87 m_nTruthCuts(truth_selection_tool.nCuts())
88 {
89 if constexpr(DetailEnabled) {
90 m_truthParticlesWithAssociatedTrack.reserve(track_to_truth_size);
91 }
92
93 };
94 void fill([[maybe_unused]] unsigned int eta_category_i,
95 [[maybe_unused]] unsigned int pdg_id_category_i,
96 [[maybe_unused]] float hit_efficiency,
97 [[maybe_unused]] float hit_purity,
98 [[maybe_unused]] float match_prob,
99 [[maybe_unused]] const xAOD::TruthParticle *best_match) {
100 BaseStat<DetailEnabled>::fill(eta_category_i,
101 pdg_id_category_i,
102 hit_efficiency,
103 hit_purity,
104 match_prob);
105 if constexpr(DetailEnabled) {
106 if (!m_truthParticlesWithAssociatedTrack.insert(best_match).second) {
107 // truth particle already had a best match
110 }
111 else {
112 ++this->m_counterPerEta[eta_category_i][kNParticleWithAssociatedTrack];
113 ++this->m_counterPerPdgId[pdg_id_category_i][kNParticleWithAssociatedTrack];
114 }
115 }
116 }
117 using TruthParticleSet = std::conditional< DetailEnabled,
118 std::unordered_set<const xAOD::TruthParticle *>,
121
126
127 unsigned int m_nTruthCuts;
128 };
130
139
143 analyseTrackTruth(const TruthParticleHitCounts &truth_particle_hit_counts,
144 const HitCountsPerTrack &track_hit_counts,
145 EventStat &event_stat) const;
146
147 void postProcessEventStat(const TruthParticleHitCounts &truth_particle_hit_counts,
148 std::size_t n_tracks,
149 EventStat &event_stat) const;
150
151 private:
153 {this, "TruthParticleHitCounts","", "Map from truth particle to hit counts." };
154
155 Gaudi::Property<std::vector<float> > m_weightsForProb
156 {this, "MatchWeights", {}, "Weights applied to the counts per measurement type for weighted sums"
157 " which are used compute the match probability." };
158 Gaudi::Property<std::vector<float> > m_weights
159 {this, "CountWeights", {}, "Weights applied to the counts per measurement type for weighted sums"
160 " which are used to compute hit efficiencies and purities." };
161
162 // Empty struct emulating the interface of a Gaudi property to replace optional properties if disabled.
163 template <class Base>
165 template <class OWNER>
166 DummyProperty( OWNER*, std::string, Base&&, std::string ) {}
167 Base value() const { return Base{}; }
168 operator Base() const { return Base{}; }
169
170 // some dummy implementations for vector properties:
171 std::size_t size() const {return 0u; }
172 bool empty() const {return true; }
173 auto operator[](std::size_t /*idx*/) { throw std::out_of_range("DummyProperty");}
174
175 // Delegate operator() to the value
176 template <class... Args>
177 decltype( std::declval<Base>()( std::declval<Args&&>()... ) ) operator()( Args&&... args ) const
178 noexcept( noexcept( std::declval<Base>()( std::declval<Args&&>()... ) ) ) {
179 return value()( std::forward<Args>( args )... );
180 }
181 };
182
183 template <class Base>
184 using Property = std::conditional< TrackFindingValidationDetailedStat, Gaudi::Property<Base>, DummyProperty<Base> >::type;
185
187 {this, "StatisticEtaBins", {-4, -2.6, -2, 0, 2., 2.6, 4}, "Gather statistics separately for these eta bins."};
189 {this, "StatisticPtBins", {1.e3,2.5e3,10e3, 100e3}, "Gather statistics separately for these pt bins."};
191 {this, "PdgIdCategorisation", false, "Categorise by pdg id."};
193 {this, "ShowRawCounts", false, "Show all counters."};
195 {this, "ShowDetailedTables", false, "Show more details; stat. uncert., RMS, entries"};
197 {this, "ComputeTrackRecoEfficiency", true, "Compute and print track reconstruction efficiency."};
198
199 ToolHandle<IAthSelectionTool> m_truthSelectionTool{this, "TruthSelectionTool","AthTruthSelectionTool", "Truth selection tool (for efficiencies and resolutions)"};
200
201
202 // helper struct for compile time optional statistics
203 template <bool IsDebug>
205 struct Empty {
206 template <typename... T_Args>
207 Empty(T_Args... ) {}
208 };
209 mutable typename std::conditional<IsDebug,
210 std::mutex,
212 mutable typename std::conditional<IsDebug,
214 Empty>::type m_measPerTruthParticleWithoutCounts ATLAS_THREAD_SAFE {20,-.5,20.-.5};
215 mutable typename std::conditional<IsDebug,
217 Empty>::type m_bestMatchProb ATLAS_THREAD_SAFE {20,0.,1.};
218 mutable typename std::conditional<IsDebug,
220 Empty>::type m_nextToBestMatchProb ATLAS_THREAD_SAFE {20,0.,1.};
221
222 template <class T_OutStream>
223 void dumpStatistics(T_OutStream &out) const;
224 void fillMeasForTruthParticleWithoutCount(double weighted_measurement_sum) const;
225 void fillTruthMatchProb(const std::array<float,2> &best_match_prob) const;
226 };
228
229 // s_NMeasurementTypes is equal to the number of UncalibMeasType
230 constexpr static unsigned int s_NMeasurementTypes = static_cast<unsigned int>(xAOD::UncalibMeasType::nTypes);
231
232 // statistics counter
242
256
257 constexpr static int s_pdgIdMax = 1000000000; // categorise all truth particles with this or larger PDG ID as "Other"
258
259
260 bool m_useAbsEtaForStat = false;
261
262 template <bool DetailEnabled>
263 struct BaseStat {
264 BaseStat() = default;
265 BaseStat([[maybe_unused]] const IAthSelectionTool &truth_selection_tool,
266 [[maybe_unused]] std::size_t per_eta_size,
267 [[maybe_unused]] std::size_t per_pdg_size)
268 : m_truthSelectionCuts(truth_selection_tool.nCuts()+1, -0.5,truth_selection_tool.nCuts()+.5)
269 {
270 if constexpr(DetailEnabled) {
271 m_counterPerEta.resize(per_eta_size);
272 m_counterPerPdgId.resize( per_pdg_size);
273 m_statPerEta.resize( per_eta_size );
274 m_statPerPdgId.resize( per_pdg_size );
275 }
276
277 }
278 void reset(const IAthSelectionTool &truth_selection_tool,
279 [[maybe_unused]] std::size_t per_eta_size,
280 [[maybe_unused]] std::size_t per_pdg_size)
281 {
282 m_truthSelectionCuts.setBinning(truth_selection_tool.nCuts()+1, -0.5,truth_selection_tool.nCuts()+.5);
283 if constexpr(DetailEnabled) {
284 m_counterPerEta.clear();
285 m_counterPerPdgId.clear();
286 m_statPerEta.clear();
287 m_statPerPdgId.clear();
288 m_counterPerEta.resize(per_eta_size);
289 m_counterPerPdgId.resize( per_pdg_size);
290 m_statPerEta.resize( per_eta_size );
291 m_statPerPdgId.resize( per_pdg_size );
292 }
293 }
294 void fill([[maybe_unused]] unsigned int eta_category_i,
295 [[maybe_unused]] unsigned int pdg_id_category_i,
296 [[maybe_unused]] float hit_efficiency,
297 [[maybe_unused]] float hit_purity,
298 [[maybe_unused]] float match_prob) {
299 if (DetailEnabled) {
300 assert( eta_category_i <m_statPerEta.size());
301 m_statPerEta[eta_category_i][kHitEfficiency].add( hit_efficiency);
302 m_statPerEta[eta_category_i][kHitPurity].add( hit_purity);
303 m_statPerEta[eta_category_i][kMatchProbability].add( match_prob);
304 assert( pdg_id_category_i <m_statPerPdgId.size());
305 m_statPerPdgId[pdg_id_category_i][kHitEfficiency].add( hit_efficiency);
306 m_statPerPdgId[pdg_id_category_i][kHitPurity].add( hit_purity);
307 m_statPerPdgId[pdg_id_category_i][kMatchProbability].add( match_prob);
308 assert( eta_category_i < m_counterPerEta.size());
309 assert( pdg_id_category_i <m_counterPerPdgId.size());
310 ++m_counterPerEta[eta_category_i][kNTotalTracks];
311 ++m_counterPerPdgId[pdg_id_category_i][kNTotalTracks];
312 }
313 }
314 void incrementTotal([[maybe_unused]] unsigned int eta_category_i,
315 [[maybe_unused]] unsigned int pdg_id_category_i) {
316 if constexpr(DetailEnabled) {
317 ++m_counterPerEta[eta_category_i][kNTotalParticles];
318 ++m_counterPerPdgId[pdg_id_category_i][kNTotalParticles];
319 }
320 }
321
323
325 const std::vector<float> &statPtBins,
326 const std::vector<float> &statEtaBins,
327 std::vector< int > &pdgId,
328 bool printDetails,
329 bool pdgIdCategorisation,
330 bool useAbsEtaForStat);
331
332
333 std::size_t perEtaSize() const {
334 if constexpr(DetailEnabled) { return m_counterPerEta.size(); }
335 else { return 0u; }
336 }
337 std::size_t perPdgIdSize() const {
338 if constexpr(DetailEnabled) { return m_counterPerPdgId.size(); }
339 else { return 0u; }
340 }
342 // per event statistics
343 struct Empty {};
344 using CounterArrayVec = std::conditional< DetailEnabled,
345 std::vector< std::array< std::size_t, kNCategorisedCounter> >,
346 Empty >::type;
347 using StatArrayVec = std::conditional< DetailEnabled,
348 std::vector< std::array<ActsUtils::Stat, kNCategorisedStat> >,
349 Empty >::type;
354
355 };
356
357
358 mutable std::mutex m_statMutex ATLAS_THREAD_SAFE;
359 mutable std::array< std::size_t, kNCounter > m_counter ATLAS_THREAD_SAFE {};
360 mutable std::vector< int > m_pdgId ATLAS_THREAD_SAFE;
362
365 void checkBinOrder( const std::vector<float> &bin_edges, const std::string &bin_label) const;
366
371 std::size_t getPtEtaStatCategory(float pt, float eta) const;
372
377 std::size_t getPtPdgIdStatCategory(float pt, int pdg_id) const;
378 void initStatTables();
379 void printStatTables() const;
380 void printCategories(const std::vector<std::string> &pt_category_labels,
381 const std::vector<std::string> &eta_category_labels,
382 std::vector<std::string> &counter_labels,
383 std::vector< std::array< ActsUtils::Stat, kNCategorisedStat> > &stat_per_category,
384 std::vector< std::array< std::size_t, kNCategorisedCounter> > &counts_per_category,
385 const std::string &top_left_label,
386 bool print_sub_categories) const;
387 void printData2D(const std::vector<std::string> &row_category_labels,
388 const std::vector<std::string> &col_category_labels,
389 const std::string &top_left_label,
390 std::vector< std::array< ActsUtils::Stat, kNCategorisedStat> > &stat_per_category,
391 std::vector< std::array< std::size_t, kNCategorisedCounter> > &counts_per_category,
392 bool rotate) const;
393
394 StatusCode checkMatchWeights();
395
396 static double weightedCountSum(const ActsTrk::HitCounterArray &counts,
397 const std::vector<float> &weights);
398 static double noiseCorrection(const ActsTrk::HitCounterArray &noise_counts,
399 const std::vector<float> &weights);
400
401 };
402
403} // namespace
404
405#endif
Scalar eta() const
pseudorapidity method
#define ATH_MSG_ERROR(x)
header file for interface of selection tools in this package
Property holding a SG store/key/clid from which a ReadHandle is made.
Property holding a SG store/key/clid from which a WriteHandle is made.
void rotate(double angler, GeoTrf::Vector2D &vector)
Container for hit counts per track Contains hit counts per associated truth particle and the total hi...
void checkBinOrder(const std::vector< float > &bin_edges, const std::string &bin_label) const
check that bins are in increasing order.
static constexpr unsigned int s_NMeasurementTypes
void printData2D(const std::vector< std::string > &row_category_labels, const std::vector< std::string > &col_category_labels, const std::string &top_left_label, std::vector< std::array< ActsUtils::Stat, kNCategorisedStat > > &stat_per_category, std::vector< std::array< std::size_t, kNCategorisedCounter > > &counts_per_category, bool rotate) const
const IAthSelectionTool & truthSelectionTool() const
std::conditional< TrackFindingValidationDetailedStat, Gaudi::Property< Base >, DummyProperty< Base > >::type Property
SG::ReadHandleKey< TruthParticleHitCounts > m_truthHitCounts
TruthMatchResult analyseTrackTruth(const TruthParticleHitCounts &truth_particle_hit_counts, const HitCountsPerTrack &track_hit_counts, EventStat &event_stat) const
DebugCounter< TrackFindingValidationDebugHists > m_debugCounter
void printCategories(const std::vector< std::string > &pt_category_labels, const std::vector< std::string > &eta_category_labels, std::vector< std::string > &counter_labels, std::vector< std::array< ActsUtils::Stat, kNCategorisedStat > > &stat_per_category, std::vector< std::array< std::size_t, kNCategorisedCounter > > &counts_per_category, const std::string &top_left_label, bool print_sub_categories) const
const TruthParticleHitCounts & getTruthParticleHitCounts(const EventContext &ctx) const
EventStatBase< TrackFindingValidationDetailedStat > EventStat
static double noiseCorrection(const ActsTrk::HitCounterArray &noise_counts, const std::vector< float > &weights)
std::size_t getPtPdgIdStatCategory(float pt, int pdg_id) const
Return the category based on the PDG ID.
static double weightedCountSum(const ActsTrk::HitCounterArray &counts, const std::vector< float > &weights)
Property< std::vector< float > > m_statPtBins
void postProcessEventStat(const TruthParticleHitCounts &truth_particle_hit_counts, std::size_t n_tracks, EventStat &event_stat) const
Gaudi::Property< std::vector< float > > m_weightsForProb
std::size_t getPtEtaStatCategory(float pt, float eta) const
Return the category based on the provided eta value.
TrackTruthMatchingBaseAlg(const std::string &name, ISvcLocator *pSvcLocator)
Property< std::vector< float > > m_statEtaBins
Gaudi::Property< std::vector< float > > m_weights
ToolHandle< IAthSelectionTool > m_truthSelectionTool
Extend Stat helper by an equidistant binned histogram.
Definition StatUtils.h:81
An algorithm that can be simultaneously executed in multiple threads.
H5Mergers.
IAthSelectionTool is a virtual baseclass for selection methods.
virtual unsigned int nCuts() const =0
return the number of cuts.
Property holding a SG store/key/clid from which a ReadHandle is made.
virtual bool isValid() override final
Can the handle be successfully dereferenced?
The AlignStoreProviderAlg loads the rigid alignment corrections and pipes them through the readout ge...
std::unordered_map< const xAOD::TruthParticle *, HitCounterArray > TruthParticleHitCounts
constexpr bool TrackFindingValidationDetailedStat
constexpr bool TrackFindingValidationDebugHists
SG::ReadCondHandle< T > makeHandle(const SG::ReadCondHandleKey< T > &key, const EventContext &ctx=Gaudi::Hive::currentContext())
TruthParticle_v1 TruthParticle
Typedef to implementation.
std::conditional< DetailEnabled, std::vector< std::array< ActsUtils::Stat, kNCategorisedStat > >, Empty >::type StatArrayVec
void reset(const IAthSelectionTool &truth_selection_tool, std::size_t per_eta_size, std::size_t per_pdg_size)
void incrementTotal(unsigned int eta_category_i, unsigned int pdg_id_category_i)
BaseStat< DetailEnabled > & operator+=(const BaseStat< DetailEnabled > &event_stat)
std::conditional< DetailEnabled, std::vector< std::array< std::size_t, kNCategorisedCounter > >, Empty >::type CounterArrayVec
void fill(unsigned int eta_category_i, unsigned int pdg_id_category_i, float hit_efficiency, float hit_purity, float match_prob)
BaseStat(const IAthSelectionTool &truth_selection_tool, std::size_t per_eta_size, std::size_t per_pdg_size)
void fillMeasForTruthParticleWithoutCount(double weighted_measurement_sum) const
std::conditional< IsDebug, std::mutex, Empty >::type m_mutex ATLAS_THREAD_SAFE
void fillTruthMatchProb(const std::array< float, 2 > &best_match_prob) const
DummyProperty(OWNER *, std::string, Base &&, std::string)
decltype(std::declval< Base >()(std::declval< Args && >()...)) operator()(Args &&... args) const noexcept(noexcept(std::declval< Base >()(std::declval< Args && >()...)))
EventStatBase(const IAthSelectionTool &truth_selection_tool, std::size_t per_eta_size, std::size_t per_pdg_size, std::size_t track_to_truth_size)
void fill(unsigned int eta_category_i, unsigned int pdg_id_category_i, float hit_efficiency, float hit_purity, float match_prob, const xAOD::TruthParticle *best_match)
std::conditional< DetailEnabled, std::unordered_set< const xAOD::TruthParticle * >, typename BaseStat< DetailEnabled >::Empty >::type TruthParticleSet
float m_hitPurity
fraction of hits originting from best match over total reco hits
float m_matchProbability
the matching probability based on weighted hit sums
float m_hitEfficiency
fraction of hits originting from best match over total best match hits
const xAOD::TruthParticle * m_truthParticle
best matching truth particle or nullptr