ATLAS Offline Software
Loading...
Searching...
No Matches
TrigTauCaloRoiUpdater.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
6
7#include "GaudiKernel/MsgStream.h"
8#include "GaudiKernel/IToolSvc.h"
9#include "GaudiKernel/StatusCode.h"
10
12#include "CxxUtils/phihelper.h"
13
14#include "TLorentzVector.h"
15
16TrigTauCaloRoiUpdater::TrigTauCaloRoiUpdater(const std::string & name, ISvcLocator* pSvcLocator)
17 : AthReentrantAlgorithm(name, pSvcLocator)
18{
19
20}
21
22
24{
25 ATH_MSG_DEBUG("Initializing " << name());
26 ATH_MSG_DEBUG("dRForCenter: " << m_dRForCenter);
27
28 ATH_MSG_DEBUG("Initialising HandleKeys");
29 ATH_CHECK(m_roIInputKey.initialize());
30 ATH_CHECK(m_clustersKey.initialize());
31 ATH_CHECK(m_roIOutputKey.initialize());
32
33 return StatusCode::SUCCESS;
34}
35
36
37StatusCode TrigTauCaloRoiUpdater::execute(const EventContext& ctx) const
38{
39 ATH_MSG_DEBUG("Running " << name());
40
41 //---------------------------------------------------------------
42 // Prepare I/O
43 //---------------------------------------------------------------
44
45 // Prepare output RoI container
46 std::unique_ptr<TrigRoiDescriptorCollection> roiCollection = std::make_unique<TrigRoiDescriptorCollection>();
48 ATH_CHECK(outputRoIHandle.record(std::move(roiCollection)));
49
50
51 // Retrieve input RoI descriptor
53 ATH_MSG_DEBUG("Size of roisHandle: " << roisHandle->size());
54 const TrigRoiDescriptor* roiDescriptor = roisHandle->at(0); // We only have one RoI in the handle
55
56
57 // Fill local variables for RoI reference position
58 float eta = roiDescriptor->eta();
59 float phi = roiDescriptor->phi();
60
61 const float dEta = (roiDescriptor->etaPlus() - roiDescriptor->etaMinus()) / 2;
62 const float dPhi = CxxUtils::deltaPhi(roiDescriptor->phiPlus(), roiDescriptor->phiMinus()) / 2;
63
64 ATH_MSG_DEBUG("RoI ID: " << roiDescriptor->roiId() << ", eta: " << eta << ", phi: " << phi);
65
66
67
68 //---------------------------------------------------------------
69 // Find detector tau axis
70 //---------------------------------------------------------------
71
72 // Retrieve Input CaloClusterContainer
74 ATH_CHECK(CCContainerHandle.isValid());
75 const xAOD::CaloClusterContainer *RoICaloClusterContainer = CCContainerHandle.get();
76
77 if(!RoICaloClusterContainer) {
78 ATH_MSG_ERROR("No CaloCluster container found");
79 return StatusCode::FAILURE;
80 }
81
82 ATH_MSG_DEBUG("Size of vector CaloCluster container is: " << RoICaloClusterContainer->size());
83
84 // We first need to get the barycenter of the LCTopo jet, including all clusters
85 TLorentzVector tau_barycenter;
86 for(const xAOD::CaloCluster* cluster : *RoICaloClusterContainer) {
87 // Skip clusters with negative energy
88 if(cluster->e() < 0) continue;
89
90 tau_barycenter += cluster->p4();
91 }
92
93 // Determine the LCTopo jet pT at the detector axis
94 TLorentzVector tau_detector_axis;
95 for(const xAOD::CaloCluster* cluster : *RoICaloClusterContainer) {
96 // Skip clusters with negative energy
97 if(cluster->e() < 0) continue;
98
99 // Skip clusters further than a maximum Delta R
100 if(tau_barycenter.DeltaR(cluster->p4()) > m_dRForCenter) continue;
101
102 tau_detector_axis += cluster->p4();
103 }
104
105
106
107 //---------------------------------------------------------------
108 // Update the RoI
109 //---------------------------------------------------------------
110
111 // Only update the roi if tau_detector_axis.Pt() > 0, i.e. if the calo cluster sum makes sense
112 if(tau_detector_axis.Pt() > 0) {
113 eta = tau_detector_axis.Eta();
114 phi = tau_detector_axis.Phi();
115 }
116
117 // Create the new RoI
118 outputRoIHandle->push_back(std::make_unique<TrigRoiDescriptor>(
119 roiDescriptor->roiWord(), roiDescriptor->l1Id(), roiDescriptor->roiId(),
120 eta, eta-dEta, eta+dEta,
122 roiDescriptor->zed(), roiDescriptor->zedMinus(), roiDescriptor->zedPlus()
123 ));
124
125
126 ATH_MSG_DEBUG("Input RoI: " << *roiDescriptor);
127 ATH_MSG_DEBUG("Output RoI: " << *outputRoIHandle->back());
128
129 return StatusCode::SUCCESS;
130}
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_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 etaMinus() const override final
gets eta at zMinus
virtual double etaPlus() const override final
gets eta at zedPlus
virtual double zed() const override final
virtual double phi() const override final
Methods to retrieve data members.
virtual double phiMinus() const override final
gets phiMinus
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 double phiPlus() const override final
gets phiPlus
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 &ctx) const override
virtual StatusCode initialize() override
SG::ReadHandleKey< TrigRoiDescriptorCollection > m_roIInputKey
SG::ReadHandleKey< xAOD::CaloClusterContainer > m_clustersKey
SG::WriteHandleKey< TrigRoiDescriptorCollection > m_roIOutputKey
TrigTauCaloRoiUpdater(const std::string &, ISvcLocator *)
Gaudi::Property< float > m_dRForCenter
T wrapToPi(T phi)
Wrap angle in radians to [-pi, pi].
Definition phihelper.h:24
T deltaPhi(T phiA, T phiB)
Return difference phiA - phiB in range [-pi, pi].
Definition phihelper.h:42
SG::ReadCondHandle< T > makeHandle(const SG::ReadCondHandleKey< T > &key, const EventContext &ctx=Gaudi::Hive::currentContext())
CaloCluster_v1 CaloCluster
Define the latest version of the calorimeter cluster class.
CaloClusterContainer_v1 CaloClusterContainer
Define the latest version of the calorimeter cluster container.
Helper for azimuthal angle calculations.