ATLAS Offline Software
Loading...
Searching...
No Matches
TrigTauRecMerged.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
5#include "TrigTauRecMerged.h"
6
7#include "GaudiKernel/SystemOfUnits.h"
12
14
17
18#include "xAODJet/Jet.h"
21
24#include "xAODTau/TauJet.h"
25
26#include "xAODTau/TauDefs.h"
27#include "xAODTau/TauTrack.h"
30
31#include "Math/Vector3D.h"
32#include "CxxUtils/phihelper.h"
33
34#include <iterator>
35#include <algorithm> //std::min
36#include <functional> //std::ref
37#include <cmath> //std::log10
38
39
40TrigTauRecMerged::TrigTauRecMerged(const std::string& name, ISvcLocator* pSvcLocator)
41 : AthReentrantAlgorithm(name, pSvcLocator)
42{
43
44}
45
46
48{
49 ATH_MSG_DEBUG("Initialize");
50
51 for(const auto& tool : m_commonTools) ATH_CHECK(tool.retrieve());
52 for(const auto& tool : m_commonToolsBeforeTF) ATH_CHECK(tool.retrieve());
53 for(const auto& tool : m_vertexFinderTools) ATH_CHECK(tool.retrieve());
54 for(const auto& tool : m_trackFinderTools) ATH_CHECK(tool.retrieve());
55 for(const auto& tool : m_vertexVarsTools) ATH_CHECK(tool.retrieve());
56 for(const auto& tool : m_idTools) ATH_CHECK(tool.retrieve());
57
58 if(!m_monTool.name().empty()) ATH_CHECK(m_monTool.retrieve());
59
60 ATH_MSG_DEBUG("Initialising handle keys");
61 ATH_CHECK(m_roiInputKey.initialize());
66
68 ATH_CHECK(m_tauJetOutputKey.initialize());
69 ATH_CHECK(m_tauTrackOutputKey.initialize());
70
71 // Hits decorations
72 if(!m_tauJetInputKey.empty() && !m_hitsInputDecorKey.empty()) {
73 // Initialize output decoration keys with the same names
74 // We manually copy the decorations because we want to make the scheduler aware of the hits dependencies
77 ATH_CHECK(m_hitsOutputDecorKey.initialize());
78
79 // Initialize input decorator keys
81 ATH_CHECK(m_hitsInputDecorKey.initialize());
82 }
83
84
85 // Propagate decorated quantities from the input to the output tau collection
86 if(!m_tauJetInputKey.empty() && !m_tauJetInputDecorKeysArray.empty()) {
87 for(auto& decor_key : m_tauJetInputDecorKeysArray) {
88 m_tauJetOutputDecorKeysArray.push_back(m_tauJetOutputKey.key() + "." + decor_key.key());
89 ATH_CHECK(m_tauJetOutputDecorKeysArray.back().initialize());
90
91 decor_key = m_tauJetInputKey.key() + "." + decor_key.key();
92 ATH_CHECK(decor_key.initialize());
93 }
94
95 } else if(m_tauJetInputKey.empty() && !m_tauJetInputDecorKeysArray.empty()) {
96 ATH_MSG_ERROR("Cannot specify input TauJet decoration keys without an input TauJet container key!");
97 return StatusCode::FAILURE;
98 }
99
100
101 // Determine the reconstruction mode
102 if(m_tauJetInputKey.empty() && !m_clustersInputKey.empty()) {
104 ATH_MSG_INFO("Running in Calo reconstruction mode (from clusters)");
105 } else if(!m_tauJetInputKey.empty()) {
107 ATH_MSG_INFO("Running in Full reconstruction mode (from input TauJets + TauTracks)");
108 } else {
109 ATH_MSG_ERROR("Invalid configuration!");
110 return StatusCode::FAILURE;
111 }
112
113
114 // Set up the monitoring accessors
115 for(const auto& [key, p] : m_monitoredIdScores) {
117 key,
118 std::make_pair(SG::ConstAccessor<float>(p.first), SG::ConstAccessor<float>(p.second))
119 );
120 }
121
122 for(const auto& [key, p] : m_monitoredHitZRegressions) {
124 key,
125 std::make_pair(SG::ConstAccessor<float>(p.first), SG::ConstAccessor<float>(p.second))
126 );
127 }
128
129
130 return StatusCode::SUCCESS;
131}
132
133
134StatusCode TrigTauRecMerged::execute(const EventContext& ctx) const
135{
136 //===============================================================================
137 // Initialize monitoring variables
138 //===============================================================================
139
140 // Common CaloOnly and Precision monitored variables:
141 auto n_taus = Monitored::Scalar<int>("NTauCandidates", 0);
142
143 auto pT = Monitored::Scalar<float>("Pt", 0);
144 auto eta = Monitored::Scalar<float>("Eta", -99.9);
145 auto phi = Monitored::Scalar<float>("Phi", -99.9);
146
147 auto etaRoI = Monitored::Scalar<float>("EtaRoI", -99.9);
148 auto phiRoI = Monitored::Scalar<float>("PhiRoI", -99.9);
149 auto dEta_RoI = Monitored::Scalar<float>("dEtaTau_RoI", -10);
150 auto dPhi_RoI = Monitored::Scalar<float>("dPhiTau_RoI", -10);
151
152 auto mEflowApprox = Monitored::Scalar<float>("mEflowApprox", -99.9);
153 auto ptRatioEflowApprox = Monitored::Scalar<float>("ptRatioEflowApprox", -99.9);
154 auto pt_jetseed_log = Monitored::Scalar<float>("pt_jetseed_log", -99.9);
155 auto etaDetectorAxis = Monitored::Scalar<float>("etaDetectorAxis", -99.9);
156 auto ptDetectorAxis = Monitored::Scalar<float>("ptDetectorAxis", -99.9);
157 auto ptDetectorAxis_log = Monitored::Scalar<float>("ptDetectorAxis_log", -99.9);
158
159 auto n_cells = Monitored::Scalar<int>("NCaloCells", 0);
160 auto EtHad = Monitored::Scalar<float>("EtHad", -10);
161 auto EtEm = Monitored::Scalar<float>("EtEm", -10);
162 auto EMFrac = Monitored::Scalar<float>("EMFrac", -10);
163 auto IsoFrac = Monitored::Scalar<float>("IsoFrac", -1);
164 auto CentFrac = Monitored::Scalar<float>("CentFrac", -10);
165
166 auto clustersMeanCenterLambda = Monitored::Scalar<float>("clustersMeanCenterLambda", 0);
167 auto clustersMeanFirstEngDens = Monitored::Scalar<float>("clustersMeanFirstEngDens", 0);
168 auto clustersMeanEMProbability = Monitored::Scalar<float>("clustersMeanEMProbability", 0);
169 auto clustersMeanSecondLambda = Monitored::Scalar<float>("clustersMeanSecondLambda", 0);
170 auto clustersMeanPresamplerFrac = Monitored::Scalar<float>("clustersMeanPresamplerFrac", 0);
171
172 auto n_clusters = Monitored::Scalar<int>("NClusters", 0);
173 std::vector<float> cluster_et_log, cluster_dEta, cluster_dPhi;
174 std::vector<float> cluster_log_SECOND_R, cluster_SECOND_LAMBDA, cluster_CENTER_LAMBDA;
175 auto mon_cluster_et_log = Monitored::Collection("cluster_et_log", cluster_et_log);
176 auto mon_cluster_dEta = Monitored::Collection("cluster_dEta", cluster_dEta);
177 auto mon_cluster_dPhi = Monitored::Collection("cluster_dPhi", cluster_dPhi);
178 auto mon_cluster_log_SECOND_R = Monitored::Collection("cluster_log_SECOND_R", cluster_log_SECOND_R);
179 auto mon_cluster_SECOND_LAMBDA = Monitored::Collection("cluster_SECOND_LAMBDA", cluster_SECOND_LAMBDA);
180 auto mon_cluster_CENTER_LAMBDA = Monitored::Collection("cluster_CENTER_LAMBDA", cluster_CENTER_LAMBDA);
181 std::vector<unsigned char> calo_errors;
182 auto mon_calo_errors = Monitored::Collection("calo_errors", calo_errors);
183
184 // Precision monitored variables
185 auto n_tracks = Monitored::Scalar<int>("NTracks", -10);
186 auto n_iso_tracks = Monitored::Scalar<int>("NIsoTracks", -10);
187
188 auto ipSigLeadTrk = Monitored::Scalar<float>("ipSigLeadTrk", -1000);
189 auto trFlightPathSig = Monitored::Scalar<float>("trFlightPathSig", -10);
190 auto massTrkSys = Monitored::Scalar<float>("massTrkSys", -10);
191 auto dRmax = Monitored::Scalar<float>("dRmax", -10);
192 auto trkAvgDist = Monitored::Scalar<float>("TrkAvgDist", -1);
193 auto innerTrkAvgDist = Monitored::Scalar<float>("innerTrkAvgDist", -1);
194 auto etovPtLead = Monitored::Scalar<float>("EtovPtLead", -10);
195 auto PSSFraction = Monitored::Scalar<float>("PSSFraction", -999.9);
196 auto EMPOverTrkSysP = Monitored::Scalar<float>("EMPOverTrkSysP", -999.9);
197 auto ChPiEMEOverCaloEME = Monitored::Scalar<float>("ChPiEMEOverCaloEME", -999.9);
198 auto vertex_x = Monitored::Scalar<float>("vertex_x", -999.9);
199 auto vertex_y = Monitored::Scalar<float>("vertex_y", -999.9);
200 auto vertex_z = Monitored::Scalar<float>("vertex_z", -999.9);
201
202 auto n_all_tracks = Monitored::Scalar<int>("NAllTracks", 0);
203 std::vector<float> track_pt_log, track_dEta, track_dPhi;
204 std::vector<float> track_d0_abs_log, track_z0sinthetaTJVA_abs_log;
205 std::vector<float> track_nPixelHitsPlusDeadSensors, track_nSCTHitsPlusDeadSensors;
206 auto mon_track_pt_log = Monitored::Collection("track_pt_log", track_pt_log);
207 auto mon_track_dEta = Monitored::Collection("track_dEta", track_dEta);
208 auto mon_track_dPhi = Monitored::Collection("track_dPhi", track_dPhi);
209 auto mon_track_d0_abs_log = Monitored::Collection("track_d0_abs_log", track_d0_abs_log);
210 auto mon_track_z0sinthetaTJVA_abs_log = Monitored::Collection("track_z0sinthetaTJVA_abs_log", track_z0sinthetaTJVA_abs_log);
211 auto mon_track_nPixelHitsPlusDeadSensors = Monitored::Collection("track_nPixelHitsPlusDeadSensors", track_nPixelHitsPlusDeadSensors);
212 auto mon_track_nSCTHitsPlusDeadSensors = Monitored::Collection("track_nSCTHitsPlusDeadSensors", track_nSCTHitsPlusDeadSensors);
213 std::vector<unsigned char> track_errors;
214 auto mon_track_errors = Monitored::Collection("track_errors", track_errors);
215
216 std::map<std::string, Monitored::Scalar<float>> monitoredIdVariables;
217 for(const auto& [key, p] : m_monitoredIdScores) {
218 monitoredIdVariables.emplace(key + "_TauJetScore_0p", Monitored::Scalar<float>(key + "_TauJetScore_0p", -1));
219 monitoredIdVariables.emplace(key + "_TauJetScoreTrans_0p", Monitored::Scalar<float>(key + "_TauJetScoreTrans_0p", -1));
220 monitoredIdVariables.emplace(key + "_TauJetScore_1p", Monitored::Scalar<float>(key + "_TauJetScore_1p", -1));
221 monitoredIdVariables.emplace(key + "_TauJetScoreTrans_1p", Monitored::Scalar<float>(key + "_TauJetScoreTrans_1p", -1));
222 monitoredIdVariables.emplace(key + "_TauJetScore_mp", Monitored::Scalar<float>(key + "_TauJetScore_mp", -1));
223 monitoredIdVariables.emplace(key + "_TauJetScoreTrans_mp", Monitored::Scalar<float>(key + "_TauJetScoreTrans_mp", -1));
224 }
225
226 // CaloHits monitored variables
227 auto n_hits = Monitored::Scalar<int>("NHits", 0);
228 std::vector<float> hit_z, hit_dPhi;
229 auto mon_hit_z = Monitored::Collection("hit_z", hit_z);
230 auto mon_hit_dPhi = Monitored::Collection("hit_dPhi", hit_dPhi);
231
232 std::map<std::string, Monitored::Scalar<float>> monitoredHitZVariables;
233 for(const auto& [key, p] : m_monitoredHitZRegressions) {
234 monitoredHitZVariables.emplace(key + "_z0", Monitored::Scalar<float>(key + "_z0", -999));
235 monitoredHitZVariables.emplace(key + "_z0_sigma", Monitored::Scalar<float>(key + "_z0_sigma", -999));
236 }
237
238 std::vector<std::reference_wrapper<Monitored::IMonitoredVariable>> monVars = {
239 std::ref(n_taus),
240 std::ref(pT), std::ref(eta), std::ref(phi),
241 std::ref(etaRoI), std::ref(phiRoI), std::ref(dEta_RoI), std::ref(dPhi_RoI),
242 std::ref(mEflowApprox), std::ref(ptRatioEflowApprox), std::ref(pt_jetseed_log),
243 std::ref(etaDetectorAxis), std::ref(ptDetectorAxis), std::ref(ptDetectorAxis_log),
244 std::ref(n_cells),
245 std::ref(EtHad), std::ref(EtEm), std::ref(EMFrac), std::ref(IsoFrac), std::ref(CentFrac),
246 std::ref(clustersMeanCenterLambda), std::ref(clustersMeanFirstEngDens), std::ref(clustersMeanEMProbability),
247 std::ref(clustersMeanSecondLambda), std::ref(clustersMeanPresamplerFrac),
248 std::ref(n_clusters), std::ref(mon_cluster_et_log), std::ref(mon_cluster_dEta), std::ref(mon_cluster_dPhi),
249 std::ref(mon_cluster_log_SECOND_R), std::ref(mon_cluster_SECOND_LAMBDA), std::ref(mon_cluster_CENTER_LAMBDA),
250 std::ref(mon_calo_errors),
251 std::ref(n_tracks), std::ref(n_iso_tracks),
252 std::ref(ipSigLeadTrk), std::ref(trFlightPathSig), std::ref(massTrkSys), std::ref(dRmax), std::ref(trkAvgDist), std::ref(innerTrkAvgDist),
253 std::ref(etovPtLead), std::ref(PSSFraction), std::ref(EMPOverTrkSysP), std::ref(ChPiEMEOverCaloEME),
254 std::ref(vertex_x), std::ref(vertex_y), std::ref(vertex_z),
255 std::ref(n_all_tracks), std::ref(mon_track_pt_log), std::ref(mon_track_dEta), std::ref(mon_track_dPhi),
256 std::ref(mon_track_d0_abs_log), std::ref(mon_track_z0sinthetaTJVA_abs_log),
257 std::ref(mon_track_nPixelHitsPlusDeadSensors), std::ref(mon_track_nSCTHitsPlusDeadSensors),
258 std::ref(mon_track_errors),
259 std::ref(n_hits), std::ref(mon_hit_z), std::ref(mon_hit_dPhi)
260 };
261 for(auto& [key, var] : monitoredIdVariables) monVars.push_back(std::ref(var));
262 for(auto& [key, var] : monitoredHitZVariables) monVars.push_back(std::ref(var));
263 auto monitorIt = Monitored::Group(m_monTool, monVars);
264
265
266 ATH_MSG_DEBUG("Executing TrigTauRecMerged");
267
268 //===============================================================================
269 // Main TauJet object and collections
270 //===============================================================================
271 xAOD::TauJet* tau = nullptr;
272
273 // We store the position of hits assigned to the tau for later monitoring
274 std::vector<ROOT::Math::XYZVector> tauHits;
275
276
277 //===============================================================================
278 // Retrieve RoI
279 //===============================================================================
281 if(!roisHandle.isValid()) {
282 ATH_MSG_ERROR("No RoIHandle found");
283 return StatusCode::FAILURE;
284 }
285
286 if(roisHandle->empty()) {
287 ATH_MSG_ERROR("Empty RoIHandle");
288 return StatusCode::FAILURE;
289 }
290
291 const TrigRoiDescriptor *roiDescriptor = roisHandle->at(0);
292 if(roiDescriptor) {
293 ATH_MSG_DEBUG("RoI: " << *roiDescriptor);
294 } else {
295 ATH_MSG_ERROR("Failed to find TrigRoiDescriptor");
296 calo_errors.push_back(NoROIDescr);
297 return StatusCode::FAILURE;
298 }
299
300
301 //===============================================================================
302 // Prepare Output TauJet and TauTrack containers
303 //===============================================================================
304
305 // Create and register the output TauJetContainer
306 std::unique_ptr<xAOD::TauJetContainer> outputContainer = std::make_unique<xAOD::TauJetContainer>();
307 std::unique_ptr<xAOD::TauJetAuxContainer> outputAuxContainer = std::make_unique<xAOD::TauJetAuxContainer>();
308 outputContainer->setStore(outputAuxContainer.get());
309
311 ATH_CHECK(outputTauHandle.record(std::move(outputContainer), std::move(outputAuxContainer)));
312
313 // Create and register the output TauTrackContainer
314 std::unique_ptr<xAOD::TauTrackContainer> outputTrackContainer = std::make_unique<xAOD::TauTrackContainer>();
315 std::unique_ptr<xAOD::TauTrackAuxContainer> outputTrackAuxContainer = std::make_unique<xAOD::TauTrackAuxContainer>();
316 outputTrackContainer->setStore(outputTrackAuxContainer.get());
317
319 ATH_CHECK(tauTrackHandle.record(std::move(outputTrackContainer), std::move(outputTrackAuxContainer)));
320
321
322 //===============================================================================
323 // Initial TauJet calo-reco / input-retrieval
324 //===============================================================================
325 // We now have two options for the TauJet reconstruction:
326 // 1) Reconstruct the TauJet from scratch, using the bare calo-clusters (1st trigger step)
327 // 2) Fetch a preceding TauJet (and [dummy] TauTracks), and use them to seed the current reconstruction
328
329 //-------------------------------------------------------------------------------
330 // Option 1: Calorimeter-only reconstruction from clusters (CaloMVA step)
331 //-------------------------------------------------------------------------------
333 // Retrieve Calocluster container
335 ATH_CHECK(CCContainerHandle.isValid());
336
337 const xAOD::CaloClusterContainer *RoICaloClusterContainer = CCContainerHandle.get();
338 if(RoICaloClusterContainer) {
339 ATH_MSG_DEBUG("CaloCluster container found of size: " << RoICaloClusterContainer->size());
340
341 // If size is zero, don't stop just continue to produce empty TauJetCollection
342 if(RoICaloClusterContainer->empty()) calo_errors.push_back(NoClustCont);
343 } else {
344 ATH_MSG_ERROR("No CaloCluster container found");
345 calo_errors.push_back(NoClustCont);
346 return StatusCode::FAILURE;
347 }
348
349 // Also create the seed-jet containers
350 std::unique_ptr<xAOD::JetContainer> jetContainer{std::make_unique<xAOD::JetContainer>()};
351 std::unique_ptr<xAOD::JetAuxContainer> jetAuxContainer{std::make_unique<xAOD::JetAuxContainer>()};
352 jetContainer->setStore(jetAuxContainer.get());
353
355 ATH_CHECK(outputTauSeedJetHandle.record(std::move(jetContainer), std::move(jetAuxContainer)));
356
357 // And create the seed-jet object itself
358 outputTauSeedJetHandle->push_back(std::make_unique<xAOD::Jet>());
359 xAOD::Jet *jet = outputTauSeedJetHandle->back();
360
361
362 // Build the jet, also keep track of the kinematics by hand
363 // Eventually, want to use FastJet here?
364 // We are using calibrated clusters, we need to keep track of this
365 jet->setConstituentsSignalState(xAOD::JetConstitScale::CalibratedJetConstituent);
366 TLorentzVector cluster_p4, barycenter;
367 for(const xAOD::CaloCluster* cluster : *RoICaloClusterContainer) {
368 ATH_MSG_DEBUG("Cluster (e, eta, phi): (" << cluster->e() << ", " << cluster->eta() << ", " << cluster->phi() << ")");
369
370 if(cluster->e() < 0) {
371 ATH_MSG_DEBUG("Negative energy cluster is rejected");
372 continue;
373 }
374
375 cluster_p4.SetPtEtaPhiE(cluster->pt(), cluster->eta(), cluster->phi(), cluster->e());
376 jet->addConstituent(cluster);
377
378 barycenter += cluster_p4;
379 }
380
381 jet->setJetP4(xAOD::JetFourMom_t(barycenter.Pt(), barycenter.Eta(), barycenter.Phi(), barycenter.M()));
382 ATH_MSG_DEBUG("Built jet with eta: " << jet->eta() << ", phi: " << jet->phi() << ", pT: " << jet->pt() << ", E: "<< jet->e() );
383
384
385 // If we're running calo-clustering, that means we just started the HLT reco,
386 // and there's no input TauJet container to this step. Create one instead!
387 outputTauHandle->push_back(std::make_unique<xAOD::TauJet>());
388 tau = outputTauHandle->back();
389 tau->setROIWord(roiDescriptor->roiWord());
390
391 // Using the new Jet collection, setup the tau candidate structure
392 tau->setJet(outputTauSeedJetHandle.ptr(), jet);
393
394 // Fix eta, phi in case the jet's energy is negative
395 if(jet->e() <= 0) {
396 ATH_MSG_DEBUG("Changing (eta, phi) back to the RoI center due to negative energy: " << jet->e());
397 tau->setP4(tau->pt(), roiDescriptor->eta(), roiDescriptor->phi(), tau->m());
398 ATH_MSG_DEBUG("Roi: " << roiDescriptor->roiId() << ", tau eta: " << tau->eta() << ", tau phi: " << tau->phi() );
399 }
400 }
401
402
403 //-------------------------------------------------------------------------------
404 // Option 2: Use input TauJet (and TauTracks) as seeds (non calo-only reco)
405 //-------------------------------------------------------------------------------
406 else if(m_reco_mode == Mode::FromTauJet) {
407 // Retrieve input TauJet container
409 const xAOD::TauJetContainer* inputTauContainer = tauInputHandle.cptr();
410 ATH_MSG_DEBUG("Input TauJet Container size: " << inputTauContainer->size());
411
412 // Copy the input TauJets to the output container
413 ATH_CHECK(deepCopy(outputTauHandle, inputTauContainer));
414
415 // Copy associated hits
416 // Store associated hits for later monitoring
417 if(!m_hitsInputDecorKey.empty()) {
418 // Initialize hits decorators
420 ATH_CHECK(hitsInputDecorHandle.isValid());
422 ATH_CHECK(hitsOutputDecorHandle.isValid());
423
424 // We use ConstAccessors instead of ReadDecorHandles to read the hit positions
425 // because hits might come from different original containers if we're including
426 // both Pixel and Strips. It is ok, because the pre-processing pipeline required to
427 // attach the hit ELs to the TauJet includes calls to HitDecoratorAlg, so the
428 // decorations will always be there, no matter how the scheduler organizes the execution.
429 static const SG::AuxElement::ConstAccessor<float> x("HitsXRelToBeamspot");
430 static const SG::AuxElement::ConstAccessor<float> y("HitsYRelToBeamspot");
431 static const SG::AuxElement::ConstAccessor<float> z("HitsZRelToBeamspot");
432
433 for(size_t i = 0; i < inputTauContainer->size(); ++i) {
434 const xAOD::TauJet* inputTau = inputTauContainer->at(i);
435 const auto& inputHits = hitsInputDecorHandle(*inputTau);
436
438 if(el.isValid()) {
439 tauHits.push_back(ROOT::Math::XYZVector(x(**el), y(**el), z(**el)));
440 }
441 }
442
443 const xAOD::TauJet* outputTau = outputTauHandle->at(i);
444 hitsOutputDecorHandle(*outputTau) = inputHits;
445 }
446 }
447
448
449 // Copy input TauJet decorations to output TauJets
450 for(size_t i = 0; i < m_tauJetInputDecorKeysArray.size(); ++i) {
452 ATH_CHECK(inputDecorHandle.isValid());
454 ATH_CHECK(outputDecorHandle.isValid());
455
456 for(size_t j = 0; j < inputTauContainer->size(); ++j) {
457 outputDecorHandle(*outputTauHandle->at(j)) = inputDecorHandle(*inputTauContainer->at(j));
458 }
459 }
460
461 // Retrieve input TauTrack container
462 if(!m_tauTrackInputKey.key().empty()) {
464 const xAOD::TauTrackContainer* inputTauTrackContainer = tauTrackInputHandle.cptr();
465 ATH_MSG_DEBUG("Tau Track Container Size " << inputTauTrackContainer->size());
466
467 // Copy the input TauTracks to the output container
468 ATH_CHECK(deepCopy(tauTrackHandle, inputTauTrackContainer));
469 }
470
471
472 // Now retrieve the main TauJet object, if available (the recently created copy of the input TauJet)
473 if(!outputTauHandle->empty()) {
474 tau = outputTauHandle->back();
475
476 // Check if the tau has a valid jetLink
477 ATH_CHECK(tau->jetLink().isValid());
478
479 // Clear all previous TauTrack links (if any), since we will run the
480 // TauTrackFinder tool (and the associated InDet helpers) in the
481 // Presel and Precision steps, refilling the tracks
482 if(!m_tauTrackInputKey.key().empty()) tau->clearTauTrackLinks();
483 }
484 }
485
486
487 //===============================================================================
488 // Get Vertex Container (optional)
489 //===============================================================================
490 const xAOD::VertexContainer* RoIVxContainer = nullptr;
491 if(!m_vertexInputKey.key().empty()){
493
494 if(!VertexContainerHandle.isValid()) {
495 ATH_MSG_DEBUG("No VertexContainer retrieved for the trigger element");
496 track_errors.push_back(NoVtxCont);
497 } else {
498 RoIVxContainer = VertexContainerHandle.get();
499 ATH_MSG_DEBUG("Size of VertexContainer: " << RoIVxContainer->size());
500 }
501 }
502
503
504 ATH_MSG_DEBUG("roiDescriptor roiWord: " << roiDescriptor->roiWord() << ", saved in TauJet: " << tau->ROIWord());
505
506
507 //===============================================================================
508 // Loop over all booked tau tools:
509 //===============================================================================
510 StatusCode processStatus = StatusCode::SUCCESS;
511
512 // Sequence: VertexFinderTools -> CommonToolsBeforeTF -> TrackFinderTools -> CommonTools -> VertexVarsTools -> IDTools
513
514 for(const auto& tool : m_vertexFinderTools) {
515 ATH_MSG_DEBUG("Starting Tool: " << tool->name());
516
517 processStatus = tool->executeVertexFinder(*tau, RoIVxContainer);
518
519 if(!processStatus.isFailure()) {
520 ATH_MSG_DEBUG(" " << tool->name() << " executed successfully");
521 } else {
522 ATH_MSG_DEBUG(" " << tool->name() << " execution failed");
523 break;
524 }
525 }
526
527 for(const auto& tool : m_commonToolsBeforeTF) {
528 if(!processStatus.isFailure()) ATH_MSG_DEBUG("Starting Tool: " << tool->name());
529 else break;
530
531 processStatus = tool->execute(*tau);
532
533 if(!processStatus.isFailure()) {
534 ATH_MSG_DEBUG(" " << tool->name() << " executed successfully");
535 } else {
536 ATH_MSG_DEBUG(" " << tool->name() << " execution failed");
537 break;
538 }
539 }
540
541 for(const auto& tool : m_trackFinderTools) {
542 if(!processStatus.isFailure()) ATH_MSG_DEBUG("Starting Tool: " << tool->name());
543 else break;
544
545 processStatus = tool->executeTrackFinder(*tau, *tauTrackHandle);
546
547 if(!processStatus.isFailure()) {
548 ATH_MSG_DEBUG(" " << tool->name() << " executed successfully");
549 } else {
550 ATH_MSG_DEBUG(" " << tool->name() << " execution failed");
551 break;
552 }
553 }
554
555 for(const auto& tool : m_commonTools) {
556 if(!processStatus.isFailure()) ATH_MSG_DEBUG("Starting Tool: " << tool->name());
557 else break;
558
559 processStatus = tool->execute(*tau);
560
561 if(!processStatus.isFailure()) {
562 ATH_MSG_DEBUG(" " << tool->name() << " executed successfully");
563 } else {
564 ATH_MSG_DEBUG(" " << tool->name() << " execution failed");
565 break;
566 }
567 }
568
569 // Dummy container passed to TauVertexVariables, not used in trigger though
570 xAOD::VertexContainer dummyVxCont;
571 for(const auto& tool : m_vertexVarsTools) {
572 if(!processStatus.isFailure()) ATH_MSG_DEBUG("Starting Tool: " << tool->name());
573 else break;
574
575 processStatus = tool->executeVertexVariables(*tau, dummyVxCont);
576
577 if(!processStatus.isFailure()) {
578 ATH_MSG_DEBUG(" " << tool->name() << " executed successfully");
579 } else {
580 ATH_MSG_DEBUG(" " << tool->name() << " execution failed");
581 break;
582 }
583 }
584
585 for(const auto& tool : m_idTools) {
586 if(!processStatus.isFailure()) ATH_MSG_DEBUG("Starting Tool: " << tool->name());
587 else break;
588
589 processStatus = tool->execute(*tau);
590
591 if(!processStatus.isFailure()) {
592 ATH_MSG_DEBUG(" " << tool->name() << " executed successfully");
593 } else {
594 ATH_MSG_DEBUG(" " << tool->name() << " execution failed");
595 break;
596 }
597 }
598
599 ATH_MSG_DEBUG("This tau has " << tau->allTracks() << " tracks linked");
600
601
602 // Cleanup in case any of the tools failed (rejected Tau)
603 if(!processStatus.isSuccess()) {
604 ATH_MSG_DEBUG("The tau object has NOT been registered in the tau container");
605
606 xAOD::TauJet* bad_tau = outputTauHandle->back();
607 ATH_MSG_DEBUG("Deleting " << bad_tau->nAllTracks() << " tracks associated with tau");
608 tauTrackHandle->erase(tauTrackHandle->end() - bad_tau->nAllTracks(), tauTrackHandle->end());
609
610 outputTauHandle->pop_back();
611
612 ATH_MSG_DEBUG("Clean up done after jet seed");
613
614 } else {
615 // Check that the seed-jet energy is positive
616 // Otherwise, try to salvage it by replacing the tau position with the RoI's
617 float fJetEnergy = (*tau->jetLink())->e();
618 ATH_MSG_DEBUG("Seed jet E: " << fJetEnergy);
619
620 if(fJetEnergy < 0.00001) {
621 ATH_MSG_DEBUG("Changing tau's (eta,phi) to RoI ones due to negative energy (PxPyPzE flips eta and phi)");
622 ATH_MSG_DEBUG("This is probably not needed anymore, method PxPyPzE has been corrected");
623 // TODO: Do we still need this??
624
625 tau->setP4(tau->pt(), roiDescriptor->eta(), roiDescriptor->phi(), tau->m());
626
627 ATH_MSG_DEBUG("Roi: " << roiDescriptor->roiId() << ", Tau eta: " << tau->eta() << ", phi: " << tau->phi() << ", pT: " << tau->pt());
628 }
629
630
631 //===============================================================================
632 // Monitor tau variables
633 //===============================================================================
634
635 pT = tau->pt()/Gaudi::Units::GeV;
636 eta = tau->eta();
637 phi = tau->phi();
638
639
640 etaRoI = roiDescriptor->eta();
641 phiRoI = roiDescriptor->phi();
642 dEta_RoI = eta - roiDescriptor->eta();
643 dPhi_RoI = phi - roiDescriptor->phi();
644 if(dPhi_RoI < -M_PI) dPhi_RoI += 2.0*M_PI;
645 if(dPhi_RoI > M_PI) dPhi_RoI -= 2.0*M_PI;
646
647
648 float pre_mEflowApprox;
649 tau->detail(xAOD::TauJetParameters::mEflowApprox, pre_mEflowApprox);
650 mEflowApprox = std::log10(std::max(pre_mEflowApprox, 140.0f));
651
652 float pre_ptRatioEflowApprox;
653 tau->detail(xAOD::TauJetParameters::ptRatioEflowApprox, pre_ptRatioEflowApprox);
654 ptRatioEflowApprox = std::min(pre_ptRatioEflowApprox, 4.0f);
655
656 pt_jetseed_log = std::log10(tau->ptJetSeed());
657 etaDetectorAxis = tau->etaDetectorAxis();
658 ptDetectorAxis = std::min(tau->ptDetectorAxis()/Gaudi::Units::GeV, 10000.0);
659 ptDetectorAxis_log = std::log10(std::min(tau->ptDetectorAxis()/Gaudi::Units::GeV, 10000.0));
660
661
664 EtHad /= Gaudi::Units::GeV;
666 EtEm /= Gaudi::Units::GeV;
667
668 float Et_raw = EtEm + EtHad;
669 if(Et_raw != 0) EMFrac = EtEm / Et_raw;
670
673
674
675 // Monitor BRT variables
676 float tmp = 0;
678 if(test) clustersMeanCenterLambda = tmp;
680 if(test) clustersMeanFirstEngDens = tmp;
682 if(test) clustersMeanEMProbability = tmp;
684 if(test) clustersMeanSecondLambda = tmp;
686 if(test) clustersMeanPresamplerFrac = tmp;
687
688
689 // Cluster variables monitoring
690 n_clusters = tau->clusters().size();
691 for(const auto& cluster : tau->clusters()) {
692 const xAOD::CaloCluster* cls = dynamic_cast<const xAOD::CaloCluster*>(cluster);
693
694 cluster_et_log.push_back(std::log10(cls->et()));
695 cluster_dEta.push_back(cls->eta() - tau->eta());
696 cluster_dPhi.push_back(cls->p4().DeltaPhi(tau->p4()));
697
698 double log_second_R = -999;
699 const bool success_SECOND_R = cls->retrieveMoment(xAOD::CaloCluster::MomentType::SECOND_R, log_second_R);
700 if(success_SECOND_R) log_second_R = std::log10(log_second_R + 0.1);
701 cluster_log_SECOND_R.push_back(log_second_R);
702
703 double second_lambda = -999;
704 const bool success_SECOND_LAMBDA = cls->retrieveMoment(xAOD::CaloCluster::MomentType::SECOND_LAMBDA, second_lambda);
705 if(success_SECOND_LAMBDA) second_lambda = std::log10(second_lambda + 0.1);
706 cluster_SECOND_LAMBDA.push_back(second_lambda);
707
708 double center_lambda = -999;
709 const bool success_CENTER_LAMBDA = cls->retrieveMoment(xAOD::CaloCluster::MomentType::CENTER_LAMBDA, center_lambda);
710 if(success_CENTER_LAMBDA) center_lambda = std::log10(center_lambda + 1e-6);
711 cluster_CENTER_LAMBDA.push_back(center_lambda);
712 }
713
714
715 // Tracks summary monitoring
716 n_tracks = tau->nTracks();
717 n_iso_tracks = tau->nTracksIsolation();
719 if(tau->nTracks() > 0) ipSigLeadTrk = std::abs(tau->track(0)->d0SigTJVA()); // TODO: Is this needed?
720 tau->detail(xAOD::TauJetParameters::trFlightPathSig, trFlightPathSig);
722 massTrkSys /= Gaudi::Units::GeV;
725 tau->detail(xAOD::TauJetParameters::innerTrkAvgDist, innerTrkAvgDist);
729 tau->detail(xAOD::TauJetParameters::ChPiEMEOverCaloEME, ChPiEMEOverCaloEME);
730
731 if(tau->vertexLink().isValid() && tau->vertex() && tau->vertex()->vertexType() != xAOD::VxType::NoVtx) {
732 vertex_x = tau->vertex()->x();
733 vertex_y = tau->vertex()->y();
734 vertex_z = tau->vertex()->z();
735 }
736
737
738 // Track variables monitoring
739 n_all_tracks = tau->allTracks().size();
740 for(const xAOD::TauTrack* track : tau->allTracks()) {
741 track_pt_log.push_back(std::log10(track->pt()));
742 track_dEta.push_back(track->eta() - tau->eta());
743 track_dPhi.push_back(track->p4().DeltaPhi(tau->p4()));
744 track_d0_abs_log.push_back(std::log10(std::abs(track->track()->d0()) + 1e-6));
745 track_z0sinthetaTJVA_abs_log.push_back(track->z0sinthetaTJVA());
746
747 uint8_t pixel_hits, pixel_dead;
748 const bool success1_pixel_hits = track->track()->summaryValue(pixel_hits, xAOD::numberOfPixelHits);
749 const bool success2_pixel_dead = track->track()->summaryValue(pixel_dead, xAOD::numberOfPixelDeadSensors);
750 float nPixelHitsPlusDeadSensor = -999;
751 if(success1_pixel_hits && success2_pixel_dead) nPixelHitsPlusDeadSensor = pixel_hits + pixel_dead;
752 track_nPixelHitsPlusDeadSensors.push_back(nPixelHitsPlusDeadSensor);
753
754 uint8_t sct_hits, sct_dead;
755 const bool success1_sct_hits = track->track()->summaryValue(sct_hits, xAOD::numberOfSCTHits);
756 const bool success2_sct_dead = track->track()->summaryValue(sct_dead, xAOD::numberOfSCTDeadSensors);
757 float nSCTHitsPlusDeadSensors = -999;
758 if(success1_sct_hits && success2_sct_dead) nSCTHitsPlusDeadSensors = sct_hits + sct_dead;
759 track_nSCTHitsPlusDeadSensors.push_back(nSCTHitsPlusDeadSensors);
760 }
761
762
763 // Hits summary monitoring
764 if(!m_hitsInputDecorKey.empty()) {
765 n_hits = tauHits.size();
766 for(const ROOT::Math::XYZVector& hit : tauHits) {
767 hit_z.push_back(hit.z());
768 hit_dPhi.push_back(CxxUtils::wrapToPi(hit.phi() - tau->phi()));
769 }
770 }
771
772
773 // TauID Score monitoring
774 for(const auto& [key, p] : m_monitoredIdScores) {
775 const auto& [score, scoreTrans] = m_monitoredInferenceAccessors.at(key);
776
777 if(score.isAvailable(*tau)) {
778 if(tau->nTracks() == 0) {
779 monitoredIdVariables.at(key + "_TauJetScore_0p") = score(*tau);
780 } else if(tau->nTracks() == 1) {
781 monitoredIdVariables.at(key + "_TauJetScore_1p") = score(*tau);
782 } else { // MP tau
783 monitoredIdVariables.at(key + "_TauJetScore_mp") = score(*tau);
784 }
785 }
786
787 if(scoreTrans.isAvailable(*tau)) {
788 if(tau->nTracks() == 0) {
789 monitoredIdVariables.at(key + "_TauJetScoreTrans_0p") = scoreTrans(*tau);
790 } else if(tau->nTracks() == 1) {
791 monitoredIdVariables.at(key + "_TauJetScoreTrans_1p") = scoreTrans(*tau);
792 } else { // MP tau
793 monitoredIdVariables.at(key + "_TauJetScoreTrans_mp") = scoreTrans(*tau);
794 }
795 }
796 }
797
798
799 // HitZ inference monitoring
800 for(const auto& [key, p] : m_monitoredHitZRegressions) {
801 const auto& [z0, sigma] = m_monitoredInferenceAccessors.at(key);
802
803 if(z0.isAvailable(*tau)) {
804 monitoredHitZVariables.at(key + "_z0") = z0(*tau);
805 }
806 if(sigma.isAvailable(*tau)) {
807 monitoredHitZVariables.at(key + "_z0_sigma") = sigma(*tau);
808 }
809 }
810
811
812 ++n_taus;
813
814
815 ATH_MSG_DEBUG("RoI: " << roiDescriptor->roiId()
816 << ", Tau pT (GeV): " << pT << ", Tau eta: " << eta << ", Tau phi: " << phi
817 << ", wrt RoI dEta: " << dEta_RoI << ", dPhi: " << dPhi_RoI);
818 }
819
820 //-------------------------------------------------------------------------------
821 // All done!
822 //-------------------------------------------------------------------------------
823
824 ATH_MSG_DEBUG("Output TauJetContainer size: " << outputTauHandle->size());
825 ATH_MSG_DEBUG("Output TauJetTrackContainer size: " << tauTrackHandle->size());
826
827 return StatusCode::SUCCESS;
828}
#define M_PI
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_DEBUG(x)
Header file to be included by clients of the Monitored infrastructure.
Handle class for reading a decoration on an object.
Handle class for adding a decoration to an object.
#define y
#define x
#define z
An algorithm that can be simultaneously executed in multiple threads.
const T * at(size_type n) const
Access an element, as an rvalue.
size_type size() const noexcept
Returns the number of elements in the collection.
bool empty() const noexcept
Returns true if the collection is empty.
Group of local monitoring quantities and retain correlation when filling histograms
Declare a monitored scalar variable.
virtual double phi() const override final
Methods to retrieve data members.
virtual double eta() const override final
SG::ConstAccessor< T, ALLOC > ConstAccessor
Definition AuxElement.h:570
Helper class to provide constant type-safe access to aux data.
Handle class for reading a decoration on an object.
virtual bool isValid() override final
Can the handle be successfully dereferenced?
const_pointer_type cptr()
Dereference the pointer.
const_pointer_type get() const
Dereference the pointer, but don't cache anything.
Handle class for adding a decoration to an object.
StatusCode record(std::unique_ptr< T > data)
Record a const object to the store.
pointer_type ptr()
Dereference the pointer.
nope - should be used for standalone also, perhaps need to protect the class def bits ifndef XAOD_ANA...
virtual unsigned int roiWord() const override final
virtual unsigned int roiId() const override final
these quantities probably don't need to be used any more
const ToolHandleArray< ITauToolBase > m_idTools
StatusCode deepCopy(SG::WriteHandle< DataVector< V > > &writeHandle, const DataVector< V > *oldContainer) const
virtual StatusCode execute(const EventContext &ctx) const override
const ToolHandleArray< ITauToolBase > m_commonTools
const ToolHandle< GenericMonitoringTool > m_monTool
SG::WriteDecorHandleKeyArray< xAOD::TauJetContainer > m_tauJetOutputDecorKeysArray
SG::WriteHandleKey< xAOD::TauJetContainer > m_tauJetOutputKey
Gaudi::Property< std::map< std::string, std::pair< std::string, std::string > > > m_monitoredHitZRegressions
SG::WriteHandleKey< xAOD::TauTrackContainer > m_tauTrackOutputKey
const ToolHandleArray< ITauToolBase > m_commonToolsBeforeTF
std::map< std::string, std::pair< SG::ConstAccessor< float >, SG::ConstAccessor< float > > > m_monitoredInferenceAccessors
SG::ReadHandleKey< xAOD::TauTrackContainer > m_tauTrackInputKey
virtual StatusCode initialize() override
SG::WriteHandleKey< xAOD::JetContainer > m_tauSeedOutputKey
SG::ReadDecorHandleKey< xAOD::TauJetContainer > m_hitsInputDecorKey
SG::WriteDecorHandleKey< xAOD::TauJetContainer > m_hitsOutputDecorKey
SG::ReadHandleKey< xAOD::VertexContainer > m_vertexInputKey
SG::ReadHandleKey< xAOD::CaloClusterContainer > m_clustersInputKey
SG::ReadDecorHandleKeyArray< xAOD::TauJetContainer > m_tauJetInputDecorKeysArray
const ToolHandleArray< ITauToolBase > m_vertexFinderTools
const ToolHandleArray< ITauToolBase > m_trackFinderTools
SG::ReadHandleKey< xAOD::TauJetContainer > m_tauJetInputKey
TrigTauRecMerged(const std::string &name, ISvcLocator *pSvcLocator)
const ToolHandleArray< ITauToolBase > m_vertexVarsTools
Gaudi::Property< std::map< std::string, std::pair< std::string, std::string > > > m_monitoredIdScores
SG::ReadHandleKey< TrigRoiDescriptorCollection > m_roiInputKey
@ SECOND_LAMBDA
Second Moment in .
@ SECOND_R
Second Moment in .
@ CENTER_LAMBDA
Shower depth at Cluster Centroid.
virtual double phi() const
The azimuthal angle ( ) of the particle.
void setROIWord(unsigned int)
void clearTauTrackLinks()
Remove all tracks from the tau.
unsigned int ROIWord() const
the ROIWord, in case TauJet is used in EF trigger
double ptDetectorAxis() const
virtual FourMom_t p4() const
The full 4-momentum of the particle.
Definition TauJet_v3.cxx:96
virtual double pt() const
The transverse momentum ( ) of the particle.
double etaDetectorAxis() const
const VertexLink_t & vertexLink() const
const JetLink_t & jetLink() const
void setJet(const xAOD::JetContainer *cont, const xAOD::Jet *jet)
bool detail(TauJetParameters::Detail detail, int &value) const
Get and set values of common details variables via enum.
const TauTrack * track(size_t i, TauJetParameters::TauTrackFlag flag=TauJetParameters::TauTrackFlag::classifiedCharged, int *container_index=0) const
Get the pointer to a given tauTrack associated with this tau /*container index needed by trackNonCons...
const Vertex * vertex() const
void setP4(double pt, double eta, double phi, double m)
Set methods for IParticle values.
std::vector< const IParticle * > clusters() const
double ptJetSeed() const
size_t nTracksIsolation() const
virtual double m() const
The invariant mass of the particle.
size_t nAllTracks() const
std::vector< const TauTrack * > allTracks() const
Get the v<const pointer> to all tracks associated with this tau, regardless of classification.
virtual double eta() const
The pseudorapidity ( ) of the particle.
size_t nTracks(TauJetParameters::TauTrackFlag flag=TauJetParameters::TauTrackFlag::classifiedCharged) const
float d0SigTJVA() const
float z() const
Returns the z position.
float y() const
Returns the y position.
VxType::VertexType vertexType() const
The type of the vertex.
float x() const
Returns the x position.
T wrapToPi(T phi)
Wrap angle in radians to [-pi, pi].
Definition phihelper.h:24
ValuesCollection< T > Collection(std::string name, const T &collection)
Declare a monitored (double-convertible) collection.
SG::ReadCondHandle< T > makeHandle(const SG::ReadCondHandleKey< T > &key, const EventContext &ctx=Gaudi::Hive::currentContext())
@ etHadAtEMScale
Get Hadronic energy at EM scale.
Definition TauDefs.h:196
@ isolFrac
Get isolation fraction.
Definition TauDefs.h:198
@ trkAvgDist
Get calibrated EM transverse energy (DEPRECATED since r19)
Definition TauDefs.h:214
@ etEMAtEMScale
Get EM energy at EM scale.
Definition TauDefs.h:194
@ centFrac
Get centrality fraction.
Definition TauDefs.h:200
@ dRmax
Get maximal dR of tracks associated to calo-seeded tau.
Definition TauDefs.h:226
@ NoVtx
Dummy vertex. TrackParticle was not used in vertex fit.
Jet_v1 Jet
Definition of the current "jet version".
CaloCluster_v1 CaloCluster
Define the latest version of the calorimeter cluster class.
VertexContainer_v1 VertexContainer
Definition of the current "Vertex container version".
TauTrack_v1 TauTrack
Definition of the current version.
Definition TauTrack.h:16
TauJet_v3 TauJet
Definition of the current "tau version".
TauTrackContainer_v1 TauTrackContainer
Definition of the current TauTrack container version.
CaloClusterContainer_v1 CaloClusterContainer
Define the latest version of the calorimeter cluster container.
TauJetContainer_v3 TauJetContainer
Definition of the current "taujet container version".
@ numberOfSCTDeadSensors
number of dead SCT sensors crossed [unit8_t].
@ numberOfSCTHits
number of hits in SCT [unit8_t].
@ numberOfPixelHits
these are the pixel hits, including the b-layer [unit8_t].
@ numberOfPixelDeadSensors
number of dead pixel sensors crossed [unit8_t].
@ CalibratedJetConstituent
Definition JetTypes.h:22
ROOT::Math::LorentzVector< ROOT::Math::PtEtaPhiM4D< double > > JetFourMom_t
Base 4 Momentum type for Jet.
Definition JetTypes.h:17
Helper for azimuthal angle calculations.