ATLAS Offline Software
Loading...
Searching...
No Matches
JetBalancePFlowJvtTool.cxx
Go to the documentation of this file.
1
2
3/*
4 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
5*/
6
7// JetBalancePFlowJvtTool.cxx
8// Implementation file for class JetBalancePFlowJvtTool
9// Author: Eimear Conroy <eimear.isobel.conroy@cern.ch>
10// Inherits from JetForwardPFlowJvtTool - Author: Anastasia Kotsokechagia <anastasia.kotsokechagia@cern.ch>
12
13// JetBalancePFlowJvtTool includes
15
16// Jet EDM
18
19// FastJet
20#include "fastjet/ClusterSequence.hh"
21#include "fastjet/ClusterSequenceArea.hh"
22#include <fastjet/AreaDefinition.hh>
23
24// Jet
26
28 // Public methods:
30
31 // Constructors
37
38 // Destructor
41 = default;
42
43 // Athena algtool's Hooks
46 {
47 ATH_MSG_INFO ("Initializing " << name() << "...");
48
49 if(m_FEKey.empty()){
50 ATH_MSG_ERROR("Flow Element container is empty");
51 return StatusCode::FAILURE;
52 }
53
55
60
61 ATH_CHECK(m_bjvtKey.initialize());
62 ATH_CHECK(m_bjvtRawKey.initialize());
63 ATH_CHECK(m_isQCDPUKey.initialize());
64 ATH_CHECK(m_isStochPUKey.initialize());
65
66 return StatusCode::SUCCESS;
67 }
68
69 StatusCode JetBalancePFlowJvtTool::decorate(const xAOD::JetContainer& jetCont) const {
70 std::vector<TVector2> pileupMomenta;
71 pileupMomenta=calculateVertexMomenta(&jetCont,m_pvind, m_vertices);
74 if(pileupMomenta.empty()) {
75 ATH_MSG_DEBUG( "pileupMomenta is empty, this can happen for events with no PU vertices."
76 <<" bJVT won't be computed for this event and will be set to 0 instead." );
77 for(const xAOD::Jet* jetC : jetCont) {
78 bjvtHandle(*jetC) = 1;
79 bjvtRawHandle(*jetC) = 0;
80 }
81 return StatusCode::SUCCESS;
82 }
83
84 for(const xAOD::Jet* jetC : jetCont) {
85 bjvtHandle(*jetC) = 1;
86 bjvtRawHandle(*jetC) = 0;
87
88 if (isCentralJet(jetC)){
89 double bjvt = getFJVT(jetC,pileupMomenta); //Same projection function as fJVT
90 if (bjvt>m_bjvtThresh) bjvtHandle(*jetC) = 0;
91 bjvtRawHandle(*jetC) = bjvt;
92 }
93 }
94 return StatusCode::SUCCESS;
95 }
96
97 //Redefined here: Using isQCDPUJet()
99 int pvind, int vertices) const {
100 std::vector<TVector2> pileupMomenta;
101 // -- Retrieve PV index if not provided by user
102 const std::size_t pv_index = (pvind==-1) ? getPV() : std::size_t(pvind);
103
105
106 for(const xAOD::Vertex* vx: *vxContHandle) {
107 if(vx->vertexType()!=xAOD::VxType::PriVtx && vx->vertexType()!=xAOD::VxType::PileUp) continue;
108 if(vx->index()==(size_t)pv_index) continue;
109
110 TString jname = m_jetsName.value();
111 jname += vx->index();
112
113 pflow::puJets vertjets = buildPFlowPUjets(*vx);
114 if( !vertjets.jetCont || !vertjets.jetAuxCont ){
115 ATH_MSG_WARNING(" Some issue appeared while building the pflow pileup jets for vertex "
116 << vx->index() << " (vxType = " << vx->vertexType()<<" )!" );
117 return pileupMomenta;
118 }
119
120 TVector2 vertex_met;
121 for( const xAOD::Jet *jet : *(vertjets.jetCont) ) {
122
123 // Remove jets which are close to hs
124 if (!m_includePV && hasCloseByHSjet(jet,pjets)) continue;
125
126 // Calculate vertex missing momentum
127 if (isQCDPUJet(jet) && getRpt(jet)> m_rptCut)
128 {
129 vertex_met += TVector2(jet->pt()*cos(jet->phi()),jet->pt()*sin(jet->phi()) ) ;
130 }
131 else{
132 vertex_met += TVector2(jet->jetP4(m_jetchargedp4).Pt()*cos(jet->jetP4(m_jetchargedp4).Phi()),
133 jet->jetP4(m_jetchargedp4).Pt()*sin(jet->jetP4(m_jetchargedp4).Phi()) );
134 }
135 }
136
137 pileupMomenta.push_back(vertex_met);
138 if(vertices!=-1 && int(vx->index())==vertices) break;
139 }
140 return pileupMomenta;
141 }
142
144 if (std::abs(jet->eta())>m_etaThresh) return false; //Must be central
145 if (jet->pt()<m_QCDPUMinPt || (m_QCDPUMaxPt>0 && jet->pt()>m_QCDPUMaxPt) ) return false;
146 //Other selections? EMPTopo studies included cuts on JVT and JVF
147 return true;
148 }
149
150 //Redefined here: Separately label QCD and stochastic PU
152 const xAOD::JetContainer *itpujets,const xAOD::JetContainer *ootpujets) {
153 //truthJets container taken to be truth HS - @to-Do double check this
154 //In-time and out-of-time PU truth jets given in AODs as separate containers
158
159 for(const xAOD::Jet *jet : *jets) {
160 bool ishs = false;
161 bool ispu = true;
162
163 bool isqcdpu = false;
164 bool isstochpu = true; //Default value to make logic work - changed below
165
166 for(const xAOD::Jet *tjet : *truthJets) {
167 if (tjet->p4().DeltaR(jet->p4())<0.3 && tjet->pt()>10e3) ishs = true;
168 if (tjet->p4().DeltaR(jet->p4())<0.6 && tjet->pt()>4e3) ispu = false;
169 }
170
171 if (ispu){
172 //Need to check both containers
173 for(const xAOD::Jet *itpujet : *itpujets){
174 if (itpujet->p4().DeltaR(jet->p4())<0.3 && itpujet->pt()>10e3) isqcdpu = true;
175 if (itpujet->p4().DeltaR(jet->p4())<0.6 && itpujet->pt()>10e3) isstochpu = false;
176 }
177 for(const xAOD::Jet *ootpujet : *ootpujets){
178 if (ootpujet->p4().DeltaR(jet->p4())<0.3 && ootpujet->pt()>10e3) isqcdpu = true;
179 if (ootpujet->p4().DeltaR(jet->p4())<0.6 && ootpujet->pt()>10e3) isstochpu = false;
180 }
181 }
182
183 else {
184 isstochpu = false; //Can't be labelled stochastic PU if jet is not determined to be PU
185 }
186
187 isHSHandle(*jet)=ishs;
188 isQCDPUHandle(*jet)=isqcdpu;
189 isStochPUHandle(*jet)=isstochpu;
190 }
191
192 return StatusCode::SUCCESS;
193 }
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
Defines enum to access jet attribute and associated particles/objects.
virtual StatusCode initialize() override
Dummy implementation of the initialisation function.
SG::WriteDecorHandleKey< xAOD::JetContainer > m_isStochPUKey
SG::WriteDecorHandleKey< xAOD::JetContainer > m_isQCDPUKey
bool isQCDPUJet(const xAOD::Jet *jet) const
Gaudi::Property< float > m_QCDPUMinPt
Gaudi::Property< float > m_QCDPUMaxPt
SG::WriteDecorHandleKey< xAOD::JetContainer > m_bjvtRawKey
virtual ~JetBalancePFlowJvtTool()
Destructor:
JetBalancePFlowJvtTool(const std::string &name)
Constructor with parameters:
StatusCode tagTruth(const xAOD::JetContainer *jets, const xAOD::JetContainer *truthJets, const xAOD::JetContainer *itpujets, const xAOD::JetContainer *ootpujets)
virtual StatusCode decorate(const xAOD::JetContainer &jetCont) const override
Decorate a jet collection without otherwise modifying it.
virtual std::vector< TVector2 > calculateVertexMomenta(const xAOD::JetContainer *jets, int pvind, int vertices) const override
SG::WriteDecorHandleKey< xAOD::JetContainer > m_bjvtKey
Gaudi::Property< double > m_etaThresh
bool hasCloseByHSjet(const xAOD::Jet *jet, const xAOD::JetContainer *pjets) const
Gaudi::Property< double > m_rptCut
virtual StatusCode initialize() override
Dummy implementation of the initialisation function.
Gaudi::Property< std::string > m_jetContainerName
SG::ReadHandleKey< xAOD::VertexContainer > m_vxContKey
SG::WriteDecorHandleKey< xAOD::JetContainer > m_isHSKey
SG::ReadHandleKey< xAOD::FlowElementContainer > m_FEKey
Gaudi::Property< std::string > m_jetchargedp4
pflow::puJets buildPFlowPUjets(const xAOD::Vertex &vx) const
Gaudi::Property< bool > m_includePV
float getFJVT(const xAOD::Jet *jet, const std::vector< TVector2 > &pileupMomenta) const
Gaudi::Property< int > m_pvind
bool isCentralJet(const xAOD::Jet *jet) const
Gaudi::Property< std::string > m_jetsName
Gaudi::Property< int > m_vertices
double getRpt(const xAOD::Jet *jet) const
JetForwardPFlowJvtTool(const std::string &name)
Constructor with parameters:
Handle class for adding a decoration to an object.
@ PileUp
Pile-up vertex.
@ PriVtx
Primary vertex.
Jet_v1 Jet
Definition of the current "jet version".
Vertex_v1 Vertex
Define the latest version of the vertex class.
JetContainer_v1 JetContainer
Definition of the current "jet container version".
std::shared_ptr< xAOD::JetAuxContainer > jetAuxCont
std::shared_ptr< xAOD::JetContainer > jetCont