ATLAS Offline Software
Loading...
Searching...
No Matches
DecoratePromptLeptonRNN.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5// Local
9
10// xAOD
15#include "GaudiKernel/ThreadLocalContext.h"
16
17// ROOT
18#include "TH1.h"
19
20//=============================================================================
21Prompt::DecoratePromptLeptonRNN::DecoratePromptLeptonRNN(const std::string& name, ISvcLocator* pSvcLocator):
22 AthAlgorithm(name, pSvcLocator),
23 m_histSvc ("THistSvc/THistSvc", name),
25{}
26
27//=============================================================================
29{
30 ATH_MSG_DEBUG("Initializing DecoratePromptLeptonRNN...");
31
32 // Initialize read/write handles
37
38 ATH_CHECK(m_eventHandleKey.initialize());
39
40 //
41 // Initialize tools and services
42 //
43 ATH_CHECK( m_histSvc.retrieve() );
44
45 m_decorHandleKeys.clear();
46 for(const std::string &label: m_toolRNN->getOutputLabels()) {
47 const std::string key = m_decorationPrefixRNN + label;
48
49 ATH_MSG_DEBUG("Add output RNN label: \"" << key << "\"");
50
51 m_decorNameMap.try_emplace (label, m_decorHandleKeys.size());
53 }
54 ATH_CHECK(m_decorHandleKeys.initialize());
55
56 ATH_MSG_DEBUG("inputContainerMuon=\"" << m_inputContainerLeptonKey << "\"");
57
58 //
59 // Instantiate Muon quality accessors
60 //
61 m_accessQuality = std::make_unique<SG::AuxElement::ConstAccessor<unsigned char> >("quality");
62
63 m_timerEvent.Reset();
64
65 ATH_MSG_DEBUG("DecoratePromptLeptonRNN initialized successfully.");
66
67 return StatusCode::SUCCESS;
68
69}
70
71//=============================================================================
73{
74 //
75 // Process current event
76 //
78
79 const EventContext& ctx = Gaudi::Hive::currentContext();
80
81 ATH_MSG_DEBUG("execute() - begin...");
82
84
85 //
86 // Retrieve object containers and event info
87 //
92
94
95 ATH_MSG_DEBUG("Size of LeptonContainer: " << h_leptons->size());
96
97 std::vector<decoratorFloatH_t> decors;
98 decors.reserve (m_decorHandleKeys.size());
100 decors.emplace_back (k, ctx);
101 }
102
103 //
104 // Find default Primary Vertex
105 //
106 const xAOD::Vertex *primaryVertex = nullptr;
107
108 for(const xAOD::Vertex *vertex: *h_vertices) {
109 if(vertex->vertexType() == xAOD::VxType::PriVtx) {
110 primaryVertex = vertex;
111 break;
112 }
113 }
114
115 //
116 // Collect tracks
117 //
118 for(const xAOD::IParticle *lepton: *h_leptons) {
119 //
120 // Select lepton track
121 //
122 const xAOD::TrackParticle *trackLep = nullptr;
123 const xAOD::Electron *elec = dynamic_cast<const xAOD::Electron*>(lepton);
124 const xAOD::Muon *muon = dynamic_cast<const xAOD::Muon *>(lepton);
125
126 if(elec) {
127 const xAOD::TrackParticle *bestmatchedGSFElTrack = elec->trackParticle(0);
128 if(bestmatchedGSFElTrack) {
129 trackLep = xAOD::EgammaHelpers::getOriginalTrackParticleFromGSF(bestmatchedGSFElTrack);
130 }
131 }
132 else if (muon) {
133 trackLep = findMuonTrack(muon);
134 }
135 else {
136 ATH_MSG_WARNING("execute - failed to find electron or muon: should never happen!");
137 }
138
139 //
140 // Find closest track jet
141 //
142 const xAOD::Jet *trackJet = findClosestTrackJet(trackLep, *h_trackJets);
143
144 if(!trackLep || !trackJet) {
145 compDummy(*lepton, decors);
146 continue;
147 }
148
149 //
150 // Select tracks within cone around lepton track.
151 //
152 std::vector<Prompt::VarHolder > select_tracks;
153
154 Prompt::VarHolder lepton_obj;
155
156 if(!prepTrackObject(lepton_obj, *trackLep, *trackLep, *trackJet, *primaryVertex, *event_handle)) {
157 continue;
158 }
159
160 //
161 // Add lepton track as one of the cone tracks.
162 //
163 select_tracks.push_back(lepton_obj);
164
165 for(const xAOD::TrackParticle *track: *h_tracks) {
166
167 if(!track) {
168 ATH_MSG_WARNING("Prompt::DecoratePromptLeptonRNN::execute - skip null track pointer - should never happen");
169 continue;
170 }
171
172 Prompt::VarHolder track_obj;
173
174 if(!prepTrackObject(track_obj, *track, *trackLep, *trackJet, *primaryVertex, *event_handle)) {
175 continue;
176 }
177
178 if(passTrack(track_obj)) {
179 select_tracks.push_back(track_obj);
180
181 ATH_MSG_DEBUG("Prompt::DecoratePromptLeptonRNN::execute - passed track pT= " << track->pt());
182 }
183 }
184
185 //
186 // Sort tracks by DR distance to lepton
187 //
188 std::sort(select_tracks.begin(), select_tracks.end(), Prompt::SortObjectByVar(Def::LepTrackDR, msg()));
189
190 //
191 // Compute RNN
192 //
193 compScore(*lepton, select_tracks, decors);
194
195 ATH_MSG_DEBUG("DecoratePromptLeptonRNN::CompScore - " << std::endl
196 << "lepton pT= " << lepton->pt()
197 << ", number of tracks: " << select_tracks.size());
198 }
199
200 return StatusCode::SUCCESS;
201}
202
203//=============================================================================
205{
206 //
207 // Finalize output
208 //
209 if(m_printTime) {
210 //
211 // Print full time stopwatch
212 //
213 ATH_MSG_INFO("DecoratePromptLeptonRNN - total time: " << PrintResetStopWatch(m_timerEvent));
214 ATH_MSG_INFO("DecoratePromptLeptonRNN - processed " << m_countEvent << " events.");
215 }
216
217 return StatusCode::SUCCESS;
218}
219
220//=============================================================================
222{
223 //
224 // Process muon - return true if all information present for RNN
225 //
226 if(muon->muonType() != xAOD::Muon::Combined || !muon->inDetTrackParticleLink().isValid()) {
227 return 0;
228 }
229
230 const unsigned char quality = (*m_accessQuality)(*muon);
231
232 ATH_MSG_DEBUG("muon pT=" << muon->pt() << " quality=" << int(quality) << " medium=" << int(xAOD::Muon::Medium));
233
234 const xAOD::TrackParticle *trackLep = *(muon->inDetTrackParticleLink());
235
236 return trackLep;
237}
238
239//=============================================================================
241 const xAOD::JetContainer &trackJets)
242{
243 //
244 // Find track jet closest to IParticle and within fixed cone
245 //
246 if(!particle) {
247 return 0;
248 }
249
250 const xAOD::Jet *trackJet = 0;
251 double currTrackJetDR = -1.0;
252
253 //
254 // Select track jet within cone around muon
255 //
256 for(const xAOD::Jet *jet: trackJets) {
257 const double dr = particle->p4().DeltaR(jet->p4());
258
259 if(currTrackJetDR < 0.0 || dr < currTrackJetDR) {
260 trackJet = jet;
261 currTrackJetDR = dr;
262 }
263 }
264
265 if(trackJet && currTrackJetDR < m_maxLepTrackJetDR) {
266 return trackJet;
267 }
268
269 return 0;
270}
271
272//=============================================================================
275 const xAOD::TrackParticle &track,
276 const xAOD::TrackParticle &lepton,
277 const xAOD::Jet &trackJet,
278 const xAOD::Vertex &priVtx,
279 const xAOD::EventInfo event
280)
281{
282 //
283 // Add xAOD::TrackParticle variables to the track object VarHolder
284 //
285 uint8_t numberOfPixelHits = 0;
286 uint8_t numberOfSCTHits = 0;
287 uint8_t numberOfPixelHoles = 0;
288 uint8_t numberOfSCTHoles = 0;
289 uint8_t numberOfPixelSharedHits = 0;
290 uint8_t numberOfSCTSharedHits = 0;
291
292 if(!(track.summaryValue(numberOfPixelHits, xAOD::numberOfPixelHits))) return false;
293 if(!(track.summaryValue(numberOfSCTHits, xAOD::numberOfSCTHits))) return false;
294 if(!(track.summaryValue(numberOfPixelHoles, xAOD::numberOfPixelHoles))) return false;
295 if(!(track.summaryValue(numberOfSCTHoles, xAOD::numberOfSCTHoles))) return false;
296 if(!(track.summaryValue(numberOfPixelSharedHits, xAOD::numberOfPixelSharedHits))) return false;
297 if(!(track.summaryValue(numberOfSCTSharedHits, xAOD::numberOfSCTSharedHits))) return false;
298
299 const uint8_t NSiHits = numberOfPixelHits + numberOfSCTHits;
300 const uint8_t NSiHoles = numberOfPixelHoles + numberOfSCTHoles;
301 const float NSiShHits = float(numberOfPixelSharedHits) + float(numberOfSCTSharedHits)/2.0;
302
303 p.addVar(Def::Pt, track.pt());
304 p.addVar(Def::AbsEta, std::fabs(track.eta()));
305 p.addVar(Def::NumberOfPIXHits, numberOfPixelHits);
306 p.addVar(Def::NumberOfSCTHits, numberOfSCTHits);
307 p.addVar(Def::NumberOfSiHits, NSiHits);
308 p.addVar(Def::NumberOfSharedSiHits, NSiShHits);
309 p.addVar(Def::NumberOfSiHoles, NSiHoles);
310 p.addVar(Def::NumberOfPixelHoles, numberOfPixelHoles);
311
312 //
313 // Add lepton - jet variables to VarHolder
314 //
315 double PtFrac = -99.;
316
317 if(track.pt() > 0.0 && trackJet.pt() > 0.0) {
318 PtFrac = track.pt() / trackJet.pt();
319 }
320
321 p.addVar(Def::TrackPtOverTrackJetPt, PtFrac);
322 p.addVar(Def::TrackJetDR, track.p4().DeltaR(trackJet.p4()));
323 p.addVar(Def::LepTrackDR, track.p4().DeltaR(lepton.p4()));
324
325 //
326 // Add Impact Parameters
327 //
328 double d0_significance = -99.;
329 double Z0Sin = 0.0;
330
331 if(track.definingParametersCovMatrixVec().size() > 0 && track.definingParametersCovMatrixVec().at(0) > 0.0) {
332 d0_significance = xAOD::TrackingHelpers::d0significance(&track,
333 event.beamPosSigmaX(),
334 event.beamPosSigmaY(),
335 event.beamPosSigmaXY());
336 }
337
338 const double deltaZ0 = track.z0() + track.vz() - priVtx.z();
339 Z0Sin = deltaZ0*std::sin(track.theta());
340
341 p.addVar(Def::Z0Sin, Z0Sin);
342 p.addVar(Def::D0Sig, d0_significance);
343
344 return true;
345}
346
347//=============================================================================
349{
350 //
351 // Select cone tracks
352 //
353 if(p.getVar(Def::LepTrackDR) < m_minTrackLeptonDR) return false;
354 if(p.getVar(Def::LepTrackDR) > m_maxTrackLeptonDR) return false;
355
356 //
357 // Kinematic track selection
358 //
359 if(p.getVar(Def::Pt) < m_minTrackpT) return false;
360 if(p.getVar(Def::AbsEta) > m_maxTrackEta) return false;
361 if(std::fabs(p.getVar(Def::Z0Sin)) > m_maxTrackZ0Sin) return false;
362
363 //
364 // Hit quality track selection
365 //
366 if(p.getVar(Def::NumberOfSiHits) < m_minTrackSiHits) return false;
367 if(p.getVar(Def::NumberOfSharedSiHits) > m_maxTrackSharedSiHits) return false;
368 if(p.getVar(Def::NumberOfSiHoles) > m_maxTrackSiHoles) return false;
369 if(p.getVar(Def::NumberOfPixelHoles) > m_maxTrackPixHoles) return false;
370
371 return true;
372}
373
374//=============================================================================
376 const std::vector<Prompt::VarHolder> &tracks,
377 std::vector<decoratorFloatH_t>& decors)
378{
379 //
380 // Call the RNN tool to get the RNN prediction for the leptons and decorate the lepton with those RNN scores.
381 //
382 ATH_MSG_DEBUG("compScore - number of tracks: " << tracks.size());
383
384 for(const Prompt::VarHolder &o: tracks) {
385 ATH_MSG_DEBUG("compScore - track: LepTrackDR = " << o.getVar(Def::LepTrackDR)
386 << ", TrackJetDR = " << o.getVar(Def::TrackJetDR)
387 << ", D0Sig = " << o.getVar(Def::D0Sig)
388 << ", Z0Sin = " << o.getVar(Def::Z0Sin)
389 << ", NumberOfPIXHits = " << o.getVar(Def::NumberOfPIXHits)
390 << ", NumberOfSCTHits = " << o.getVar(Def::NumberOfSCTHits)
391 << ", PtFrac = " << o.getVar(Def::TrackPtOverTrackJetPt) );
392 }
393
394 const std::map<std::string, double> results = m_toolRNN->computeRNNOutput(tracks);
395
396 for(const std::pair<const std::string, double>& v: results) {
397 //
398 // Decorate muon
399 //
400
401 ATH_MSG_DEBUG("DecoratePromptLeptonRNN compScore - " << v.first << " = " << v.second );
402
403 auto dit = m_decorNameMap.find (v.first);
404
405 if(dit != m_decorNameMap.end()) {
406 decors.at(dit->second)(particle) = v.second;
407 }
408 else {
409 ATH_MSG_WARNING("CompScore - unknown output label=\"" << v.first << "\"");
410 }
411
412 if(m_debug) {
413 std::map<std::string, TH1*>::iterator hit = m_hists.find(v.first);
414
415 if(hit == m_hists.end()) {
416 TH1* h = 0;
417
418 StatusCode hist_status = makeHist(h, v.first, 100, 0.0, 1.0);
419 if (hist_status != StatusCode::SUCCESS){
420 ATH_MSG_WARNING("DecoratePromptLeptonRNN compScore - failed to make hist");
421 }
422
423 hit = m_hists.insert(std::map<std::string, TH1*>::value_type(v.first, h)).first;
424 }
425
426 if(hit->second) {
427 hit->second->Fill(v.second);
428 }
429 }
430 }
431
432 return true;
433}
434
435//=============================================================================
437 std::vector<decoratorFloatH_t>& decors) const
438{
439 //
440 // Fill dummy values for RNN outputs
441 //
442 for (auto& d : decors) {
443 d(particle) = -1.0;
444 }
445
446 return true;
447}
448
449//=============================================================================
450StatusCode Prompt::DecoratePromptLeptonRNN::makeHist(TH1 *&h, const std::string &key, int nbin, double xmin, double xmax)
451{
452 //
453 // Initiliase histogram pointer. If configured to run in validation mode, then create and register histogram
454 //
455 h = 0;
456
457 if(m_outputStream.empty() || key.empty()) {
458 return StatusCode::SUCCESS;
459 }
460
461 const std::string hname = name() + "_" + key;
462 const std::string hist_key = "/"+m_outputStream+"/"+hname;
463
464 h = new TH1D(hname.c_str(), hname.c_str(), nbin, xmin, xmax);
465 h->SetDirectory(0);
466
467 return m_histSvc->regHist(hist_key, h);
468}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_INFO(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
Some common helper functions used by decoration handles.
AthAlgorithm(const std::string &name, ISvcLocator *pSvcLocator)
Constructor with parameters:
MsgStream & msg() const
Header file for AthHistogramAlgorithm.
virtual StatusCode finalize() override
SG::ReadHandleKey< xAOD::JetContainer > m_inputContainerTrackJetKey
Gaudi::Property< double > m_maxTrackEta
const xAOD::TrackParticle * findMuonTrack(const xAOD::Muon *muon)
bool compDummy(const xAOD::IParticle &particle, std::vector< decoratorFloatH_t > &decors) const
Gaudi::Property< unsigned > m_minTrackSiHits
std::map< std::string, TH1 * > m_hists
DecoratePromptLeptonRNN(const std::string &name, ISvcLocator *pSvcLocator)
Gaudi::Property< unsigned > m_maxTrackSiHoles
StatusCode makeHist(TH1 *&h, const std::string &key, int nbin, double xmin, double xmax)
std::unordered_map< std::string, size_t > m_decorNameMap
Gaudi::Property< double > m_maxTrackLeptonDR
const xAOD::Jet * findClosestTrackJet(const xAOD::TrackParticle *particle, const xAOD::JetContainer &trackJets)
SG::ReadHandleKey< xAOD::IParticleContainer > m_inputContainerLeptonKey
Gaudi::Property< double > m_maxLepTrackJetDR
std::unique_ptr< SG::AuxElement::ConstAccessor< unsigned char > > m_accessQuality
Gaudi::Property< double > m_maxTrackZ0Sin
virtual StatusCode initialize() override
SG::WriteDecorHandleKeyArray< xAOD::IParticleContainer > m_decorHandleKeys
SG::ReadHandleKey< xAOD::VertexContainer > m_inputContainerPrimaryVerticesKey
Gaudi::Property< std::string > m_decorationPrefixRNN
Gaudi::Property< std::string > m_outputStream
SG::ReadHandleKey< xAOD::TrackParticleContainer > m_inputContainerTrackKey
bool compScore(const xAOD::IParticle &particle, const std::vector< Prompt::VarHolder > &tracks, std::vector< decoratorFloatH_t > &decors)
virtual StatusCode execute() override
Gaudi::Property< double > m_maxTrackSharedSiHits
Gaudi::Property< unsigned > m_maxTrackPixHoles
SG::ReadHandleKey< xAOD::EventInfo > m_eventHandleKey
Gaudi::Property< double > m_minTrackLeptonDR
bool prepTrackObject(Prompt::VarHolder &p, const xAOD::TrackParticle &track, const xAOD::TrackParticle &lepton, const xAOD::Jet &trackJet, const xAOD::Vertex &priVtx, const xAOD::EventInfo event)
Gaudi::Property< double > m_minTrackpT
Property holding a SG store/key/clid/attr name from which a WriteDecorHandle is made.
const xAOD::TrackParticle * trackParticle(size_t index=0) const
Pointer to the xAOD::TrackParticle/s that match the electron candidate.
Class providing the definition of the 4-vector interface.
virtual FourMom_t p4() const
The full 4-momentum of the particle.
Definition Jet_v1.cxx:71
virtual double pt() const
The transverse momentum ( ) of the particle.
Definition Jet_v1.cxx:44
virtual FourMom_t p4() const override final
The full 4-momentum of the particle.
float z() const
Returns the z position.
std::string label(const std::string &format, int i)
Definition label.h:19
double xmax
Definition listroot.cxx:61
double xmin
Definition listroot.cxx:60
@ NumberOfPixelHoles
Definition VarHolder.h:53
@ NumberOfPIXHits
Definition VarHolder.h:48
@ NumberOfSCTHits
Definition VarHolder.h:49
@ NumberOfSharedSiHits
Definition VarHolder.h:51
@ TrackPtOverTrackJetPt
Definition VarHolder.h:55
@ NumberOfSiHoles
Definition VarHolder.h:52
std::string PrintResetStopWatch(TStopwatch &watch)
std::string makeContDecorKey(const std::string &cont, const std::string &decor)
Make a StoreGate key from container and decoration name.
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.
const xAOD::TrackParticle * getOriginalTrackParticleFromGSF(const xAOD::TrackParticle *trkPar)
Helper function for getting the "Original" Track Particle (i.e before GSF) via the GSF Track Particle...
double d0significance(const xAOD::TrackParticle *tp, double d0_uncert_beam_spot_2)
@ PriVtx
Primary vertex.
Jet_v1 Jet
Definition of the current "jet version".
EventInfo_v1 EventInfo
Definition of the latest event info version.
TrackParticle_v1 TrackParticle
Reference the current persistent version:
Vertex_v1 Vertex
Define the latest version of the vertex class.
Muon_v1 Muon
Reference the current persistent version:
JetContainer_v1 JetContainer
Definition of the current "jet container version".
@ numberOfPixelHoles
number of pixel layers on track with absence of hits [unit8_t].
@ numberOfSCTHits
number of hits in SCT [unit8_t].
@ numberOfPixelHits
these are the pixel hits, including the b-layer [unit8_t].
@ numberOfPixelSharedHits
number of Pixel all-layer hits shared by several tracks [unit8_t].
@ numberOfSCTSharedHits
number of SCT hits shared by several tracks [unit8_t].
@ numberOfSCTHoles
number of SCT holes [unit8_t].
Electron_v1 Electron
Definition of the current "egamma version".