ATLAS Offline Software
TrigTauCaloRoiUpdater.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2022 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 
16 TrigTauCaloRoiUpdater::TrigTauCaloRoiUpdater(const std::string & name, ISvcLocator* pSvcLocator) :
17  AthReentrantAlgorithm(name, pSvcLocator) {}
18 
20 
21  ATH_MSG_DEBUG( "Initializing " << name() );
22  ATH_MSG_DEBUG( "z0HalfWidth " << m_z0HalfWidth );
23  if(m_z0HalfWidth <= 0.) {
24  ATH_MSG_DEBUG( "z0HalfWidth <= 0: will use the original RoIInput z0HalfWidth" );
25  }
26  ATH_MSG_DEBUG( "etaHalfWidth " << m_etaHalfWidth );
27  ATH_MSG_DEBUG( "phiHalfWidth " << m_phiHalfWidth );
28  ATH_MSG_DEBUG( "dRForCenter " << m_dRForCenter );
29 
30  ATH_MSG_DEBUG( "Initialising HandleKeys" );
34 
35  return StatusCode::SUCCESS;
36 }
37 
38 
39 
40 StatusCode TrigTauCaloRoiUpdater::execute(const EventContext& ctx) const {
41 
42  ATH_MSG_DEBUG( "Running "<< name() <<" ... " );
43 
44  // Prepare Outputs
45  std::unique_ptr< TrigRoiDescriptorCollection > roICollection( new TrigRoiDescriptorCollection() );
46 
47  // Retrieve Input CaloClusterContainer
49  CHECK( CCContainerHandle.isValid() );
50  const xAOD::CaloClusterContainer *RoICaloClusterContainer = CCContainerHandle.get();
51 
52  if(RoICaloClusterContainer != nullptr) {
53  ATH_MSG_DEBUG( "Size of vector CaloCluster container is " << RoICaloClusterContainer->size());
54  if(RoICaloClusterContainer->empty()) {
55  ATH_MSG_DEBUG( "CaloCluster container is empty");
56  }
57  }else {
58  ATH_MSG_ERROR( "no CaloCluster container found " );
59  return StatusCode::FAILURE;
60  }
61 
62  //get RoI descriptor
64  ATH_MSG_DEBUG("Size of roisHandle: "<<roisHandle->size());
65  const TrigRoiDescriptor *roiDescriptor = roisHandle->at(0);
66 
67  // fill local variables for RoI reference position
68  float eta = roiDescriptor->eta();
69  float phi = roiDescriptor->phi();
70  float dEta = m_etaHalfWidth;
71  float dPhi = m_phiHalfWidth;
72 
73  float zed = roiDescriptor->zed();
74  float zedPlus = roiDescriptor->zedPlus();
75  float zedMinus = roiDescriptor->zedMinus();
76  if(m_z0HalfWidth > 0.) {
77  zedPlus = zed + m_z0HalfWidth;
78  zedMinus = zed - m_z0HalfWidth;
79  }
80 
81  ATH_MSG_DEBUG( "; RoI ID = " << roiDescriptor->roiId()
82  << ": Eta = " << eta
83  << ", Phi = " << phi );
84 
85  // Make a minimal effort to speed things up ;)
86  TLorentzVector myCluster;
87  TLorentzVector TauBarycenter(0., 0., 0., 0.);
88 
90  for (clusterIt=RoICaloClusterContainer->begin(); clusterIt != RoICaloClusterContainer->end(); ++clusterIt) {
91  if((*clusterIt)->e() < 0)
92  continue;
93 
94  myCluster.SetPtEtaPhiE((*clusterIt)->pt(), (*clusterIt)->eta(), (*clusterIt)->phi(), (*clusterIt)->e());
95  TauBarycenter += myCluster;
96  }
97 
98  // Determine the LC tau pT at detector axis
99  TLorentzVector TauDetectorAxis(0.,0.,0.,0.);
100  for (clusterIt=RoICaloClusterContainer->begin(); clusterIt != RoICaloClusterContainer->end(); ++clusterIt) {
101 
102  if((*clusterIt)->e() < 0)
103  continue;
104 
105  myCluster.SetPtEtaPhiE((*clusterIt)->pt(), (*clusterIt)->eta(), (*clusterIt)->phi(), (*clusterIt)->e());
106  if(TauBarycenter.DeltaR(myCluster) > m_dRForCenter)
107  continue;
108 
109  TauDetectorAxis += myCluster;
110  } // end loop on clusters
111 
112  //Only update the roi if TauDetectorAxis.Pt() is larger than zero, in other words, if the calo clusters sum makes sense
113  if(TauDetectorAxis.Eta()!=roiDescriptor->eta() && TauDetectorAxis.Pt()>0.) eta = TauDetectorAxis.Eta();
114  if(TauDetectorAxis.Phi()!=roiDescriptor->phi() && TauDetectorAxis.Pt()>0.) phi = TauDetectorAxis.Phi();
115 
116  // Prepare the new RoI
117  TrigRoiDescriptor *outRoi = new TrigRoiDescriptor(roiDescriptor->roiWord(), roiDescriptor->l1Id(), roiDescriptor->roiId(),
118  eta, eta-dEta, eta+dEta,
120  zed, zedMinus, zedPlus);
121 
122  ATH_MSG_DEBUG("Input RoI " << *roiDescriptor);
123  ATH_MSG_DEBUG("Output RoI " << *outRoi);
124 
125  roICollection->push_back(outRoi);
126 
127  // Save Outputs
128  ATH_MSG_DEBUG( "Saving RoIs to be used as input to Fast Tracking -- TO BE CHANGED -- ::: " << m_roIOutputKey.key() );
130  CHECK( outputRoiHandle.record( std::move( roICollection ) ) );
131 
132  return StatusCode::SUCCESS;
133 }
TrigTauCaloRoiUpdater::m_dRForCenter
Gaudi::Property< float > m_dRForCenter
Definition: TrigTauCaloRoiUpdater.h:26
DataModel_detail::const_iterator
Const iterator class for DataVector/DataList.
Definition: DVLIterator.h:82
TrigTauCaloRoiUpdater::m_roIOutputKey
SG::WriteHandleKey< TrigRoiDescriptorCollection > m_roIOutputKey
Definition: TrigTauCaloRoiUpdater.h:33
phi
Scalar phi() const
phi method
Definition: AmgMatrixBasePlugin.h:64
RoiDescriptor::zedMinus
virtual double zedMinus() const override final
z at the most backward end of the RoI
Definition: RoiDescriptor.h:113
eta
Scalar eta() const
pseudorapidity method
Definition: AmgMatrixBasePlugin.h:79
SG::ReadHandle
Definition: StoreGate/StoreGate/ReadHandle.h:70
CxxUtils::wrapToPi
T wrapToPi(T phi)
Wrap angle in radians to [-pi, pi].
Definition: phihelper.h:24
TrigRoiDescriptor::roiWord
virtual unsigned int roiWord() const override final
Definition: TrigRoiDescriptor.h:135
TrigTauCaloRoiUpdater::m_phiHalfWidth
Gaudi::Property< float > m_phiHalfWidth
Definition: TrigTauCaloRoiUpdater.h:29
TrigRoiDescriptor
nope - should be used for standalone also, perhaps need to protect the class def bits #ifndef XAOD_AN...
Definition: TrigRoiDescriptor.h:56
SG::VarHandleKey::key
const std::string & key() const
Return the StoreGate ID for the referenced object.
Definition: AthToolSupport/AsgDataHandles/Root/VarHandleKey.cxx:141
AthReentrantAlgorithm
An algorithm that can be simultaneously executed in multiple threads.
Definition: AthReentrantAlgorithm.h:83
SG::makeHandle
SG::ReadCondHandle< T > makeHandle(const SG::ReadCondHandleKey< T > &key, const EventContext &ctx=Gaudi::Hive::currentContext())
Definition: ReadCondHandle.h:270
TrigTauCaloRoiUpdater::m_etaHalfWidth
Gaudi::Property< float > m_etaHalfWidth
Definition: TrigTauCaloRoiUpdater.h:28
ATH_MSG_ERROR
#define ATH_MSG_ERROR(x)
Definition: AthMsgStreamMacros.h:33
EL::StatusCode
::StatusCode StatusCode
StatusCode definition for legacy code.
Definition: PhysicsAnalysis/D3PDTools/EventLoop/EventLoop/StatusCode.h:22
SG::ReadHandle::get
const_pointer_type get() const
Dereference the pointer, but don't cache anything.
ATH_MSG_DEBUG
#define ATH_MSG_DEBUG(x)
Definition: AthMsgStreamMacros.h:29
TauGNNUtils::Variables::Track::dPhi
bool dPhi(const xAOD::TauJet &tau, const xAOD::TauTrack &track, double &out)
Definition: TauGNNUtils.cxx:530
TrigTauCaloRoiUpdater.h
TrigTauCaloRoiUpdater::TrigTauCaloRoiUpdater
TrigTauCaloRoiUpdater(const std::string &, ISvcLocator *)
Definition: TrigTauCaloRoiUpdater.cxx:16
TrigRoiDescriptorCollection
Athena::TPCnvVers::Current Athena::TPCnvVers::Old TrigRoiDescriptorCollection
Definition: TrigSteeringEventTPCnv.cxx:78
TrigRoiDescriptor::l1Id
virtual unsigned int l1Id() const override final
Definition: TrigRoiDescriptor.h:134
CHECK
#define CHECK(...)
Evaluate an expression and check for errors.
Definition: Control/AthenaKernel/AthenaKernel/errorcheck.h:422
SG::VarHandleKey::initialize
StatusCode initialize(bool used=true)
If this object is used as a property, then this should be called during the initialize phase.
Definition: AthToolSupport/AsgDataHandles/Root/VarHandleKey.cxx:103
DataVector
Derived DataVector<T>.
Definition: DataVector.h:581
SG::ReadHandle::isValid
virtual bool isValid() override final
Can the handle be successfully dereferenced?
RoiDescriptor::zed
virtual double zed() const override final
Definition: RoiDescriptor.h:102
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:195
phihelper.h
Helper for azimuthal angle calculations.
DataVector::push_back
value_type push_back(value_type pElem)
Add an element to the end of the collection.
TrigRoiDescriptor::roiId
virtual unsigned int roiId() const override final
these quantities probably don't need to be used any more
Definition: TrigRoiDescriptor.h:133
DataVector::end
const_iterator end() const noexcept
Return a const_iterator pointing past the end of the collection.
SG::WriteHandle< TrigRoiDescriptorCollection >
TrigTauCaloRoiUpdater::m_roIInputKey
SG::ReadHandleKey< TrigRoiDescriptorCollection > m_roIInputKey
Definition: TrigTauCaloRoiUpdater.h:31
TrigTauCaloRoiUpdater::m_z0HalfWidth
Gaudi::Property< float > m_z0HalfWidth
Definition: TrigTauCaloRoiUpdater.h:27
SG::WriteHandle::record
StatusCode record(std::unique_ptr< T > data)
Record a const object to the store.
RoiDescriptor::phi
virtual double phi() const override final
Methods to retrieve data members.
Definition: RoiDescriptor.h:100
TrigTauCaloRoiUpdater::initialize
virtual StatusCode initialize() override
Definition: TrigTauCaloRoiUpdater.cxx:19
RoiDescriptor::eta
virtual double eta() const override final
Definition: RoiDescriptor.h:101
TrigRoiDescriptor
Athena::TPCnvVers::Current TrigRoiDescriptor
Definition: TrigSteeringEventTPCnv.cxx:68
TrigRoiDescriptor.h
DataVector::at
const T * at(size_type n) const
Access an element, as an rvalue.
TrigTauCaloRoiUpdater::execute
virtual StatusCode execute(const EventContext &ctx) const override
Definition: TrigTauCaloRoiUpdater.cxx:40
TauGNNUtils::Variables::Track::dEta
bool dEta(const xAOD::TauJet &tau, const xAOD::TauTrack &track, double &out)
Definition: TauGNNUtils.cxx:525
DataVector::size
size_type size() const noexcept
Returns the number of elements in the collection.
DataVector::empty
bool empty() const noexcept
Returns true if the collection is empty.
DataVector::begin
const_iterator begin() const noexcept
Return a const_iterator pointing at the beginning of the collection.
TrigTauCaloRoiUpdater::m_clustersKey
SG::ReadHandleKey< xAOD::CaloClusterContainer > m_clustersKey
Definition: TrigTauCaloRoiUpdater.h:32
RoiDescriptor::zedPlus
virtual double zedPlus() const override final
z at the most forward end of the RoI
Definition: RoiDescriptor.h:112