ATLAS Offline Software
Loading...
Searching...
No Matches
TrigTauTrackRoiUpdater.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
5#include <cmath>
6
8
9#include "GaudiKernel/IToolSvc.h"
10#include "GaudiKernel/StatusCode.h"
11
13#include "CxxUtils/phihelper.h"
14
16
17TrigTauTrackRoiUpdater::TrigTauTrackRoiUpdater(const std::string& name, ISvcLocator* pSvcLocator)
18 : AthReentrantAlgorithm(name, pSvcLocator)
19{
20
21}
22
23
28
29
31{
32 ATH_MSG_DEBUG("Initializing " << name() );
33 ATH_MSG_DEBUG("z0HalfWidth: " << m_z0HalfWidth );
34 ATH_MSG_DEBUG("etaHalfWidth: " << m_etaHalfWidth );
35 ATH_MSG_DEBUG("phiHalfWidth: " << m_phiHalfWidth );
36 ATH_MSG_DEBUG("nHitPix: " << m_nHitPix );
37 ATH_MSG_DEBUG("nSiHoles: " << m_nSiHoles );
38
39 if(m_z0HalfWidth < 0 || m_etaHalfWidth < 0 || m_phiHalfWidth < 0) {
40 ATH_MSG_ERROR("Incorrect parameters");
41 return StatusCode::FAILURE;
42 }
43
44 ATH_MSG_DEBUG("Initialising HandleKeys");
45 ATH_CHECK(m_roIInputKey.initialize());
46 ATH_CHECK(m_tracksKey.initialize());
47 ATH_CHECK(m_roIOutputKey.initialize());
48
49 return StatusCode::SUCCESS;
50}
51
52
53StatusCode TrigTauTrackRoiUpdater::execute(const EventContext& ctx) const
54{
55 ATH_MSG_DEBUG("Running " << name());
56
57 //---------------------------------------------------------------
58 // Prepare I/O
59 //---------------------------------------------------------------
60
61 // Prepare output RoI container
62 std::unique_ptr<TrigRoiDescriptorCollection> roiCollection = std::make_unique<TrigRoiDescriptorCollection>();
64 ATH_CHECK(outputRoIHandle.record(std::move(roiCollection)));
65
66
67 // Retrieve Input TrackCollection
69 ATH_CHECK(TrackCollectionHandle.isValid());
70 const xAOD::TrackParticleContainer* foundTracks = TrackCollectionHandle.get();
71
72 if(!foundTracks) {
73 ATH_MSG_ERROR("No track container found, the Track RoI updater should not be scheduled");
74 return StatusCode::FAILURE;
75 }
76 ATH_MSG_DEBUG("Found " << foundTracks->size() << " FTF tracks, updating the RoI");
77
78
79 // Retrieve input RoI descriptor
81 ATH_MSG_DEBUG("Size of roisHandle: " << roisHandle->size());
82 const TrigRoiDescriptor* roiDescriptor = roisHandle->at(0); // We only have one RoI in the handle
83
84
85 // Fill local variables for RoI reference position
86 float eta = roiDescriptor->eta();
87 float phi = roiDescriptor->phi();
88
89 float zed = roiDescriptor->zed();
90 float zedMinus = roiDescriptor->zedMinus();
91 float zedPlus = roiDescriptor->zedPlus();
92
93
94
95 //---------------------------------------------------------------
96 // Find leading track
97 //---------------------------------------------------------------
98
99 const xAOD::TrackParticle* leadTrack = nullptr;
100 float trkPtMax = 0;
101
102 // Use the highest-pt track satisfying quality cuts
103 // If no track is found, the input ROI is used
104 for(const xAOD::TrackParticle* track : *foundTracks) {
105
106 float trackPt = track->pt();
107 if(trackPt > trkPtMax) {
108 uint8_t nPix{}, nPixHoles{}, nSCTHoles{}, summaryVal{};
109 nPix = track->summaryValue(summaryVal, xAOD::numberOfPixelHits) ? summaryVal : 0;
110 if(nPix < m_nHitPix) {
111 ATH_MSG_DEBUG("Track rejected because nHitPix " << static_cast<int>(nPix) << " < " << m_nHitPix);
112 continue;
113 }
114
115 nPixHoles = track->summaryValue(summaryVal, xAOD::numberOfPixelHoles) ? summaryVal : 0;
116 nSCTHoles = track->summaryValue(summaryVal, xAOD::numberOfSCTHoles) ? summaryVal : 0;
117 if((nPixHoles + nSCTHoles) > m_nSiHoles) {
118 ATH_MSG_DEBUG("Track rejected because nSiHoles " << static_cast<int>(nPixHoles + nSCTHoles) << " > " << m_nSiHoles);
119 continue;
120 }
121
122 leadTrack = track;
123 trkPtMax = trackPt;
124 ATH_MSG_VERBOSE("pTmax = " << trkPtMax);
125 }
126 }
127
128
129
130 //---------------------------------------------------------------
131 // Update the RoI
132 //---------------------------------------------------------------
133
134 // If a leading track is found, update all the ROI position and size parameters;
135 // else, only update the eta/phi width (etaMinus, etaPlus, phiMinus, phiPlus)
136 if(leadTrack) {
137 zed = leadTrack->z0() + leadTrack->vz();
138 ATH_MSG_DEBUG("Track z0 " << leadTrack->z0() <<" vz: " << leadTrack->vz() << " zed:" << zed);
139 zedMinus = zed - m_z0HalfWidth;
140 zedPlus = zed + m_z0HalfWidth;
141 eta = leadTrack->eta();
142 phi = leadTrack->phi0();
143 }
144
145 float etaMinus = eta - m_etaHalfWidth;
146 float etaPlus = eta + m_etaHalfWidth;
147 float phiMinus = CxxUtils::wrapToPi(phi - m_phiHalfWidth);
148 float phiPlus = CxxUtils::wrapToPi(phi + m_phiHalfWidth);
149
150
151 // Create the new RoI
152 outputRoIHandle->push_back(std::make_unique<TrigRoiDescriptor>(
153 roiDescriptor->roiWord(), roiDescriptor->l1Id(), roiDescriptor->roiId(),
154 eta, etaMinus, etaPlus,
155 phi, phiMinus, phiPlus,
156 zed, zedMinus, zedPlus
157 ));
158
159
160 ATH_MSG_DEBUG("Input RoI: " << *roiDescriptor);
161 ATH_MSG_DEBUG("Output RoI: " << *outputRoIHandle->back());
162
163 return StatusCode::SUCCESS;
164}
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_VERBOSE(x)
#define ATH_MSG_DEBUG(x)
An algorithm that can be simultaneously executed in multiple threads.
size_type size() const noexcept
Returns the number of elements in the collection.
virtual double zed() const override final
virtual double phi() const override final
Methods to retrieve data members.
virtual double zedPlus() const override final
z at the most forward end of the RoI
virtual double zedMinus() const override final
z at the most backward end of the RoI
virtual double eta() const override final
virtual bool isValid() override final
Can the handle be successfully dereferenced?
const_pointer_type get() const
Dereference the pointer, but don't cache anything.
StatusCode record(std::unique_ptr< T > data)
Record a const object to the store.
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
virtual unsigned int l1Id() const override final
virtual StatusCode execute(const EventContext &) const override
Gaudi::Property< int > m_nSiHoles
SG::ReadHandleKey< TrigRoiDescriptorCollection > m_roIInputKey
SG::ReadHandleKey< xAOD::TrackParticleContainer > m_tracksKey
Gaudi::Property< float > m_phiHalfWidth
SG::WriteHandleKey< TrigRoiDescriptorCollection > m_roIOutputKey
Gaudi::Property< int > m_nHitPix
Gaudi::Property< float > m_z0HalfWidth
TrigTauTrackRoiUpdater(const std::string &, ISvcLocator *)
Gaudi::Property< float > m_etaHalfWidth
virtual StatusCode initialize() override
float z0() const
Returns the parameter.
float vz() const
The z origin for the parameters.
virtual double eta() const override final
The pseudorapidity ( ) of the particle.
float phi0() const
Returns the parameter, which has range to .
T wrapToPi(T phi)
Wrap angle in radians to [-pi, pi].
Definition phihelper.h:24
SG::ReadCondHandle< T > makeHandle(const SG::ReadCondHandleKey< T > &key, const EventContext &ctx=Gaudi::Hive::currentContext())
TrackParticle_v1 TrackParticle
Reference the current persistent version:
TrackParticleContainer_v1 TrackParticleContainer
Definition of the current "TrackParticle container version".
@ numberOfPixelHoles
number of pixel layers on track with absence of hits [unit8_t].
@ numberOfPixelHits
these are the pixel hits, including the b-layer [unit8_t].
@ numberOfSCTHoles
number of SCT holes [unit8_t].
Helper for azimuthal angle calculations.