ATLAS Offline Software
Loading...
Searching...
No Matches
HIJetConstituentModifierTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3*/
4
12
14
18
20
21 ATH_MSG_INFO("Initializing HIJetConstituentModifierTool w/ Cluster Key " << m_clusterKey.key() );
22 ATH_CHECK( m_clusterKey.initialize() );
23 return StatusCode::SUCCESS;
24
25}
26
28
29 float E_min=m_subtractorTool->minEnergyForMoments();
30 const xAOD::JetConstituentVector constituents = jet.getConstituents();
31 std::vector<size_t> cluster_indices;
32 cluster_indices.reserve(constituents.size());
33
34 for (const auto *constituent : constituents)
35 {
36 cluster_indices.push_back(constituent->rawConstituent()->index());
37 }
38
41 constituentAcc( "constituentLinks" );
43 constituentWeightAcc( "constituentWeights" );
44
45 if( constituentAcc.isAvailable(jet) ) constituentAcc( jet ).resize(0);
46 if( constituentWeightAcc.isAvailable(jet) ) constituentWeightAcc( jet ).resize(0);
47
48 //save unsubtracted kinematics as moment if they don't exist already...
49 xAOD::IParticle::FourMom_t unsubtractedP4;
50 unsubtractedP4 = jet.p4();
51 jet.setJetP4(HIJetRec::unsubtractedJetState(),jet.jetP4());
52
53 xAOD::IParticle::FourMom_t subtractedP4;
54 xAOD::IParticle::FourMom_t subtractedOriginCorrP4;
55 xAOD::JetFourMom_t subtractedJet4vec;
56 xAOD::JetFourMom_t subtractedOriginCorrJet4vec;
57 //Cluster container access (from the e-gamma + cluster subtracted deep copy)
58 SG::ReadHandle<xAOD::CaloClusterContainer> readHandleSubtractedClusters ( m_clusterKey );
59
60 const xAOD::CaloClusterContainer* ccl=readHandleSubtractedClusters.get();
61
62 for(auto index : cluster_indices)
63 {
64 const xAOD::CaloCluster* cl= ccl->at(index);
65
66 jet.addConstituent(cl);
67 xAOD::IParticle::FourMom_t subtractedClusterP4;
69 subtractedP4+=subtractedClusterP4;
70 ATH_MSG_DEBUG("Subracted Cluster #: " << cl->index() << " :: E: " << subtractedClusterP4.E() << " :: Eta: " << subtractedClusterP4.Eta() << " :: Phi: " << subtractedClusterP4.Phi());
71
73 xAOD::IParticle::FourMom_t subtractedOriginCorrClusterP4;
75 subtractedOriginCorrP4+=subtractedOriginCorrClusterP4;
76 ATH_MSG_DEBUG("Subrtracted OC Cluster #: " << cl->index() << " :: E: " << subtractedOriginCorrClusterP4.E() << " :: Eta: " << subtractedOriginCorrClusterP4.Eta() << " :: Phi: " << subtractedOriginCorrClusterP4.Phi());
77 }
78
79 }
80
81 //Check for subtracted Kinematics
82 if(subtractedP4.E()/std::cosh(subtractedP4.Eta()) < E_min)
83 {
84 subtractedP4=unsubtractedP4;
85 subtractedP4*=1e-7;//ghost scale
86 }
87 //Check for subtracted + Origin Correction Kinematics
88 if(subtractedOriginCorrP4.E()/std::cosh(subtractedOriginCorrP4.Eta()) < E_min && m_originCorrection )
89 {
90 subtractedOriginCorrP4=unsubtractedP4;
91 subtractedOriginCorrP4*=1e-7;//ghost scale
92 }
93
94 subtractedJet4vec.SetCoordinates(subtractedP4.Pt(),subtractedP4.Eta(),
95 subtractedP4.Phi(),subtractedP4.M());
96 jet.setJetP4(HIJetRec::subtractedJetState(), subtractedJet4vec);
97
98 if( m_originCorrection ) {
99 subtractedOriginCorrJet4vec.SetCoordinates(subtractedOriginCorrP4.Pt(),subtractedOriginCorrP4.Eta(),
100 subtractedOriginCorrP4.Phi(),subtractedOriginCorrP4.M());
101 jet.setJetP4(HIJetRec::subtractedOriginCorrectedJetState(), subtractedOriginCorrJet4vec);
102 jet.setJetP4(subtractedOriginCorrJet4vec);
103 jet.setConstituentsSignalState(HIJetRec::subtractedOriginCorrectedConstitState());
104
105 }
106 else {
107 jet.setJetP4(subtractedJet4vec);
108 jet.setConstituentsSignalState(HIJetRec::subtractedConstitState());
109 }
110
111 return 0;
112}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_INFO(x)
#define ATH_MSG_DEBUG(x)
This file defines helper classes to deal with jet constituents.
Handle class for reading from StoreGate.
const T * at(size_type n) const
Access an element, as an rvalue.
SG::ReadHandleKey< xAOD::CaloClusterContainer > m_clusterKey
Name of input cluster container.
HIJetConstituentModifierTool(const std::string &myname)
Gaudi::Property< bool > m_originCorrection
|brief boolean switch to drive the JetScale settings after constituents are added
virtual int modifyJet(xAOD::Jet &jet) const override
Modify a single jet. This is obsolete and set to be removed.
virtual StatusCode initialize() override
Dummy implementation of the initialisation function.
ToolHandle< IHISubtractorTool > m_subtractorTool
handle to IHISubtractorTool that determines the subtracted kinematics for each constituent
JetModifierBase(const std::string &myname)
Ctor.
SG::Accessor< T, ALLOC > Accessor
Definition AuxElement.h:572
bool isAvailable(const ELT &e) const
Test to see if this variable exists in the store.
const_pointer_type get() const
Dereference the pointer, but don't cache anything.
TLorentzVector FourMom_t
Definition of the 4-momentum type.
A vector of jet constituents at the scale used during jet finding.
size_t size() const
number of constituents
constexpr const char * unsubtractedJetState()
constexpr xAOD::JetConstitScale subtractedOriginCorrectedConstitState()
constexpr const char * subtractedJetState()
TLorentzVector getClusterP4(const xAOD::CaloCluster *cl, xAOD::CaloCluster::State s)
: getClusterP4 Added during Rel22 migration.
constexpr const char * subtractedOriginCorrectedJetState()
constexpr xAOD::CaloCluster::State subtractedClusterState()
constexpr xAOD::JetConstitScale subtractedConstitState()
constexpr xAOD::CaloCluster::State subtractedOriginCorrectedClusterState()
Definition index.py:1
Jet_v1 Jet
Definition of the current "jet version".
CaloCluster_v1 CaloCluster
Define the latest version of the calorimeter cluster class.
CaloClusterContainer_v1 CaloClusterContainer
Define the latest version of the calorimeter cluster container.
ROOT::Math::LorentzVector< ROOT::Math::PtEtaPhiM4D< double > > JetFourMom_t
Base 4 Momentum type for Jet.
Definition JetTypes.h:17