ATLAS Offline Software
Loading...
Searching...
No Matches
dRMatchingTool.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 INDETPHYSVALMONITORING_DRMATCHINGTOOL_H
6#define INDETPHYSVALMONITORING_DRMATCHINGTOOL_H 1
7
14
15// STL include(s)
16#include <atomic>
17#include <cmath> /* std::fabs, sqrt, pow */
18#include <algorithm> /* std::all_of (C++11), std::sort, std::lower_bound, std::upper_bound, std::set_intersection */
19#include <iterator> /* std::back_inserter */
20#include <mutex>
21
22// ATLAS software include(s)
27#include "AsgTools/AsgTool.h"
28
29// xAOD include(s)
34
35// Local include(s)
37
85 public virtual ::IAsgSelectionTool,
86 public asg::AsgTool {
88public:
90 dRMatchingTool(const std::string& name);
91
93 virtual
95
97 virtual StatusCode initialize() override;
98 virtual StatusCode finalize() override;
99 virtual const asg::AcceptInfo& getAcceptInfo( ) const override;
100 virtual asg::AcceptData accept(const xAOD::IParticle* p) const override;
101
103 // For matching single track to set of truth particles.
104 virtual asg::AcceptData
106 const xAOD::TruthParticleContainer* truthParticles,
107 bool (* truthSelectionTool)(const xAOD::TruthParticle*) = nullptr) const;
108
111 const xAOD::TruthParticleContainer* truthParticles,
112 bool (* truthSelectionTool)(const xAOD::TruthParticle*) = nullptr) const;
113
114 // For matching single truth particle to set of tracks.
115 virtual asg::AcceptData
117 const xAOD::TrackParticleContainer* trackParticles,
118 bool (* trackSelectionTool)(const xAOD::TrackParticle*) = nullptr) const;
119
122 const xAOD::TrackParticleContainer* trackParticles,
123 bool (* trackSelectionTool)(const xAOD::TrackParticle*) = nullptr) const;
124
125 // Accessor for the minimal dR from the latest matching of either type.
126 float dRmin() const;
127protected:
128 // Struct for event cache
129 struct CacheEntry {
130 EventContext::ContextEvt_t m_evt{EventContext::INVALID_CONTEXT_EVT};
131
132 // The minimal dR for latest successful matching of either type.
133 float m_dRmin = 9999.;
134
135 // Copy of the xAOD::TruthParticleContainer used to perform the caching. A
136 // copy, rather than a pointer or reference, is necessary to check when a new
137 // container is encountered (i.e. when a new event is processed).
139 // Copy of the xAOD::TrackParticleContainer used to perform the caching.
141
142 // Cached vectors of the xAOD::TruthParticles in m_baseTruthContainer, sorted
143 // in ascending order accoring to their pT, eta, and phi, resp.
144 std::vector< const xAOD::TruthParticle* > m_truthParticlesSortedPt;
145 std::vector< const xAOD::TruthParticle* > m_truthParticlesSortedEta;
146 std::vector< const xAOD::TruthParticle* > m_truthParticlesSortedPhi;
147 // Cached vectors of the xAOD::TrackParticles in m_baseTrackContainer, sorted
148 // in ascending order accoring to their pT, eta, and phi, resp.
149 std::vector< const xAOD::TrackParticle* > m_trackParticlesSortedPt;
150 std::vector< const xAOD::TrackParticle* > m_trackParticlesSortedEta;
151 std::vector< const xAOD::TrackParticle* > m_trackParticlesSortedPhi;
152
153 void check(EventContext::ContextEvt_t evt) {
154 // Check event number and reset variables if a new event
155 if (m_evt==evt) return;
156
157 m_dRmin = 9999.;
164 m_evt = evt;
165 }
166 };
167
169 // Cache track particles.
171 CacheEntry* ent,
172 bool (* trackSelectionTool)(const xAOD::TrackParticle*) = nullptr) const;
173 // Cache truth particles.
175 CacheEntry* ent,
176 bool (* truthSelectionTool)(const xAOD::TruthParticle*) = nullptr) const;
177
178 // Clear cached track particle arrays.
179 inline void
181 ent->m_trackParticlesSortedPt.clear();
182 ent->m_trackParticlesSortedEta.clear();
183 ent->m_trackParticlesSortedPhi.clear();
184 return;
185 }
186
187 // Clear cached truth particle arrays.
188 inline void
190 ent->m_truthParticlesSortedPt.clear();
191 ent->m_truthParticlesSortedEta.clear();
192 ent->m_truthParticlesSortedPhi.clear();
193 return;
194 }
195
196 // Sort and cache all particles of type U in container into the three vectors,
197 // subject to their passing the selectionTool.
198 template<class T, class U>
199 void sortVectors(const T* container,
200 std::vector< const U* >& vec_pt,
201 std::vector< const U* >& vec_eta,
202 std::vector< const U* >& vec_phi,
203 bool (* selectionTool)(const U*)) const;
204
205 // Match the particle p by dR and/or pT to the particles in the three cached
206 // vectors.
207 template<class U, class V>
208 bool sortedMatch(const U* p,
209 std::vector< const V* >& vec_pt,
210 std::vector< const V* >& vec_eta,
211 std::vector< const V* >& vec_phi,
212 float& dRmin) const;
213private:
215 // AcceptInfo object.
217 // Vector of stored cut names and descriptions.
218 std::vector<std::pair<std::string, std::string> > m_cuts;
219
220 // Counters.
221 mutable std::atomic<ULong64_t> m_numProcessed{0}; // !< a counter of the number of tracks proccessed
222 mutable std::atomic<ULong64_t> m_numPassed{0}; // !< a counter of the number of tracks that passed all cuts
223 mutable std::vector<ULong64_t> m_numPassedCuts ATLAS_THREAD_SAFE; // !< tracks the number of tracks that passed each cut family. Guarded by m_mutex
224 mutable std::mutex m_mutex; // !< To guard m_numPassedCuts and m_cache
225
227 // The maximal dR-distance allowed for successful matching.
228 FloatProperty m_dRmax{this, "dRmax", -1.};
229 // The maximal pT-resolution allowed for successful matching.
230 FloatProperty m_pTResMax{this, "pTResMax", -1.};
231
232 // Event cache
233 mutable SG::SlotSpecificObj<CacheEntry> m_cache ATLAS_THREAD_SAFE; // Guarded by m_mutex
234};
235
236#endif // > !INDETPHYSVALMONITORING_DRMATCHINGTOOL_H
Maintain a set of objects, one per slot.
Define macros for attributes used to control the static checker.
Maintain a set of objects, one per slot.
Base class for the dual-use tool implementation classes.
Definition AsgTool.h:47
dRMatchingTool(const std::string &name)
Constructor(s).
virtual StatusCode finalize() override
FloatProperty m_dRmax
Cut vales.
ASG_TOOL_CLASS1(dRMatchingTool, IAsgSelectionTool)
std::atomic< ULong64_t > m_numProcessed
virtual const asg::AcceptInfo & getAcceptInfo() const override
Declare the interface ID for this pure-virtual interface class to the Athena framework.
void checkCacheTruthParticles(const xAOD::TruthParticleContainer *truthParticles, CacheEntry *ent, bool(*truthSelectionTool)(const xAOD::TruthParticle *)=nullptr) const
virtual asg::AcceptData accept(const xAOD::IParticle *p) const override
The main accept method: the actual cuts are applied here.
std::atomic< ULong64_t > m_numPassed
void clearTrackParticles(CacheEntry *ent) const
bool sortedMatch(const U *p, std::vector< const V * > &vec_pt, std::vector< const V * > &vec_eta, std::vector< const V * > &vec_phi, float &dRmin) const
std::vector< std::pair< std::string, std::string > > m_cuts
virtual StatusCode initialize() override
SelectionTool method(s).
std::mutex m_mutex
virtual asg::AcceptData acceptLegacy(const xAOD::TrackParticle *p, const xAOD::TruthParticleContainer *truthParticles, bool(*truthSelectionTool)(const xAOD::TruthParticle *)=nullptr) const
dR-matching specific method(s).
float dRmin() const
void checkCacheTrackParticles(const xAOD::TrackParticleContainer *trackParticles, CacheEntry *ent, bool(*trackSelectionTool)(const xAOD::TrackParticle *)=nullptr) const
Internal method(s).
std::vector< ULong64_t > m_numPassedCuts ATLAS_THREAD_SAFE
virtual ~dRMatchingTool()
Destructor.
void clearTruthParticles(CacheEntry *ent) const
asg::AcceptInfo m_accept
Data member(s).
FloatProperty m_pTResMax
void sortVectors(const T *container, std::vector< const U * > &vec_pt, std::vector< const U * > &vec_eta, std::vector< const U * > &vec_phi, bool(*selectionTool)(const U *)) const
Class providing the definition of the 4-vector interface.
TrackParticle_v1 TrackParticle
Reference the current persistent version:
TruthParticle_v1 TruthParticle
Typedef to implementation.
TrackParticleContainer_v1 TrackParticleContainer
Definition of the current "TrackParticle container version".
TruthParticleContainer_v1 TruthParticleContainer
Declare the latest version of the truth particle container.
std::vector< const xAOD::TrackParticle * > m_trackParticlesSortedPt
std::vector< const xAOD::TrackParticle * > m_trackParticlesSortedEta
xAOD::TruthParticleContainer m_baseTruthContainer
std::vector< const xAOD::TruthParticle * > m_truthParticlesSortedEta
void check(EventContext::ContextEvt_t evt)
std::vector< const xAOD::TrackParticle * > m_trackParticlesSortedPhi
EventContext::ContextEvt_t m_evt
std::vector< const xAOD::TruthParticle * > m_truthParticlesSortedPt
std::vector< const xAOD::TruthParticle * > m_truthParticlesSortedPhi
xAOD::TrackParticleContainer m_baseTrackContainer