ATLAS Offline Software
Loading...
Searching...
No Matches
JetForwardPFlowJvtTool.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// JetForwardPFlowJvtTool.cxx
8// Implementation file for class JetForwardPFlowJvtTool
9// Author: Anastasia Kotsokechagia <anastasia.kotsokechagia@cern.ch>
11
12// JetForwardPFlowJvtTool includes
14
15// Jet EDM
17
18// FastJet
19#include "fastjet/ClusterSequence.hh"
20#include "fastjet/ClusterSequenceArea.hh"
21#include <fastjet/AreaDefinition.hh>
22
23// Jet
25#include "JetRec/JetClusterer.h"
26
28 // Public methods:
30
31 // Constructors
34 AsgTool(name) {
35 }
36
37 // Destructor
40 = default;
41
42 // Athena algtool's Hooks
45 {
46 ATH_MSG_INFO ("Initializing " << name() << "...");
47 if (m_tightOP) m_fjvtThresh = 0.53;
48 else m_fjvtThresh = 0.72;
49
50 ATH_CHECK( m_tvaKey.initialize() );
51
52 if(m_jetContainerName.empty()){
53 ATH_MSG_ERROR("JetForwardPFlowJvtTool needs to have its input jet container configured!");
54 return StatusCode::FAILURE;
55 }
56
57 //FlowElement jets instead of PFO jets
58 if(!m_FEKey.empty()){
59 if(!m_orFEKey.key().empty()){
61 }
62 }
63 else{ //PFO reconstruction
64 if(!m_orKey.key().empty()){
65 m_orKey = m_jetContainerName + "." + m_orKey.key();
66 }
67 }
68 ATH_CHECK(m_PFOKey.initialize( m_FEKey.empty() ));
69 ATH_CHECK(m_orKey.initialize( m_FEKey.empty() && !m_orKey.key().empty() ));
70 ATH_CHECK(m_FEKey.initialize( !m_FEKey.empty() ));
71 ATH_CHECK(m_orFEKey.initialize( !m_FEKey.empty() && !m_orFEKey.key().empty() ));
72 ATH_CHECK(m_eventInfoKey.initialize());
73
78 m_passJvtKey = m_jetContainerName + "." + m_passJvtKey.key(); //nnjvt pass
79
80 ATH_CHECK(m_fjvtKey.initialize());
81 ATH_CHECK(m_fjvtRawKey.initialize());
82 ATH_CHECK(m_isHSKey.initialize());
83 ATH_CHECK(m_isPUKey.initialize());
84 ATH_CHECK(m_passJvtKey.initialize());
85
86 ATH_CHECK(m_vxContKey.initialize());
87
88#ifndef XAOD_STANDALONE
91 }
92#endif
93
94 return StatusCode::SUCCESS;
95 }
96
97 StatusCode JetForwardPFlowJvtTool::decorate(const xAOD::JetContainer& jetCont) const {
98 std::vector<TVector2> pileupMomenta;
99
100 pileupMomenta=calculateVertexMomenta(&jetCont,m_pvind, m_vertices);
101
104 if(pileupMomenta.empty()) {
105 ATH_MSG_DEBUG( "pileupMomenta is empty, this can happen for events with no PU vertices."
106 <<" fJVT won't be computed for this event and will be set to 0 instead." );
107 for(const xAOD::Jet* jetF : jetCont) {
108 fjvtHandle(*jetF) = 1;
109 fjvtRawHandle(*jetF) = 0;
110 }
111 return StatusCode::SUCCESS;
112 }
113
114 for(const xAOD::Jet* jetF : jetCont) {
115 fjvtHandle(*jetF) = 1;
116 fjvtRawHandle(*jetF) = 0;
117
118 if (isForwardJet(jetF)){
119 double fjvt = getFJVT(jetF,pileupMomenta);
120 if (fjvt>m_fjvtThresh) fjvtHandle(*jetF) = 0;
121 fjvtRawHandle(*jetF) = fjvt;
122 }
123 }
124 return StatusCode::SUCCESS;
125 }
126
127 float JetForwardPFlowJvtTool::getFJVT(const xAOD::Jet *jet, const std::vector<TVector2>& pileupMomenta) const {
128 TVector2 fjet(jet->pt()*cos(jet->phi()),jet->pt()*sin(jet->phi()));
129 double fjvt = 0;
130 for (const TVector2& pu : pileupMomenta) {
131 double projection = pu*fjet/fjet.Mod();
132 if (projection<fjvt) fjvt = projection;
133 }
134 return -1*fjvt/fjet.Mod();
135 }
136
138 int pvind, int vertices) const {
139 std::vector<TVector2> pileupMomenta;
140 // -- Retrieve PV index if not provided by user
141 const std::size_t pv_index = (pvind==-1) ? getPV() : std::size_t(pvind);
142
144
145 for(const xAOD::Vertex* vx: *vxContHandle) {
146 if(vx->vertexType()!=xAOD::VxType::PriVtx && vx->vertexType()!=xAOD::VxType::PileUp) continue;
147 if(vx->index()==(size_t)pv_index) continue;
148
149 TString jname = m_jetsName.value();
150 jname += vx->index();
151
152 pflow::puJets vertjets = buildPFlowPUjets(*vx);
153 if( !vertjets.jetCont || !vertjets.jetAuxCont ){
154 ATH_MSG_WARNING(" Some issue appeared while building the pflow pileup jets for vertex "
155 << vx->index() << " (vxType = " << vx->vertexType()<<" )!" );
156 return pileupMomenta;
157 }
158
159 TVector2 vertex_met;
160 for( const xAOD::Jet *jet : *(vertjets.jetCont) ) {
161
162 // Remove jets which are close to hs
163 if (!m_includePV && hasCloseByHSjet(jet,pjets)) continue;
164
165 // Calculate vertex missing momentum
167 {
168 vertex_met += TVector2(jet->pt()*cos(jet->phi()),jet->pt()*sin(jet->phi()) ) ;
169 }
170 else{
171 vertex_met += TVector2(jet->jetP4(m_jetchargedp4).Pt()*cos(jet->jetP4(m_jetchargedp4).Phi()),
172 jet->jetP4(m_jetchargedp4).Pt()*sin(jet->jetP4(m_jetchargedp4).Phi()) );
173 }
174 }
175
176 pileupMomenta.push_back(vertex_met);
177 if(vertices!=-1 && int(vx->index())==vertices) break;
178 }
179 return pileupMomenta;
180 }
181
183 for (const xAOD::Jet* pjet : *pjets) {
184 char jet_nnjvtpass=false;
186 jet_nnjvtpass = passJvtHandle(*pjet);
187 if (pjet->p4().DeltaR(jet->p4())<0.3 && jet_nnjvtpass && isCentralJet(pjet) ) return true;
188 }
189 return false;
190 }
191
193 pflow::puJets pu_jets;
194 const std::size_t pv_index = (m_pvind==-1) ? getPV() : std::size_t (m_pvind);
195
196 std::vector<fastjet::PseudoJet> input_pfo;
197 std::set<int> charged_pfo;
198
200
201 if (!tvaHandle.isValid()){
202 ATH_MSG_ERROR("Could not retrieve the TrackVertexAssociation: "
203 << m_tvaKey.key());
204 return pu_jets;
205 }
206
207 if(!m_FEKey.empty()){
209
210 for(const xAOD::FlowElement* fe : *FlowElementHandle){
211 if (!m_orFEKey.key().empty()){
213 if (!orHandle(*fe)) continue;
214 }
215 if (fe->isCharged()) {
216 const xAOD::TrackParticle* track = dynamic_cast<const xAOD::TrackParticle*>(fe->chargedObject(0));
217
218 if (vx.index()==pv_index && std::abs((vx.z()-track->z0())*sin(track->theta()))>m_dzCut)
219 continue;
220 if (vx.index()!=pv_index
221 && (!tvaHandle->associatedVertex(track)
222 || vx.index()!=tvaHandle->associatedVertex(track)->index())
223 ) continue;
224 input_pfo.push_back(feToPseudoJet(fe, CP::charged, &vx) );
225 charged_pfo.insert(fe->index());
226 }
227 else if (std::abs(fe->eta())<m_neutMaxRap && !fe->isCharged() && fe->e()>0)
228 {
229 input_pfo.push_back(feToPseudoJet(fe, CP::neutral, &vx) );
230 }
231 }
232 }
233 else{
235
236 for(const xAOD::PFO* pfo : *PFOHandle){
237 if (!m_orKey.key().empty()){
239 if (!orHandle(*pfo)) continue;
240 }
241 if (pfo->isCharged()) {
242 if (vx.index()==pv_index && std::abs((vx.z()-pfo->track(0)->z0())*sin(pfo->track(0)->theta()))>m_dzCut)
243 continue;
244 if (vx.index()!=pv_index
245 && (!tvaHandle->associatedVertex(pfo->track(0))
246 || vx.index()!=tvaHandle->associatedVertex(pfo->track(0))->index())
247 ) continue;
248 input_pfo.push_back(pfoToPseudoJet(pfo, CP::charged, &vx) );
249 charged_pfo.insert(pfo->index());
250 }
251 else if (std::abs(pfo->eta())<m_neutMaxRap && !pfo->isCharged() && pfo->eEM()>0)
252 {
253 input_pfo.push_back(pfoToPseudoJet(pfo, CP::neutral, &vx) );
254 }
255 }
256 }
257
258 std::shared_ptr<xAOD::JetContainer> vertjets = std::make_shared<xAOD::JetContainer>();
259 std::shared_ptr<xAOD::JetAuxContainer> vertjetsAux = std::make_shared<xAOD::JetAuxContainer>();
260
261 vertjets->setStore(vertjetsAux.get());
262 TString newname = m_jetsName.value();
263 newname += vx.index();
264
265 // Use run/event number as random number seeds.
266 std::vector<int> seeds;
267 auto evtInfoHandle = SG::makeHandle(m_eventInfoKey);
268 if (!evtInfoHandle.isValid())
269 {
270 ATH_MSG_ERROR("Unable to retrieve event info");
271 return pu_jets;
272 }
273 JetClustererHelper::seedsFromEventInfo(evtInfoHandle.cptr(), seeds);
274
275 fastjet::JetDefinition jet_def(fastjet::antikt_algorithm,0.4);
276 fastjet::AreaDefinition area_def(fastjet::active_area_explicit_ghosts,
277 fastjet::GhostedAreaSpec(fastjet::SelectorAbsRapMax(m_maxRap)));
278 area_def = area_def.with_fixed_seed(seeds);
279 fastjet::ClusterSequenceArea clust_pfo(input_pfo,jet_def,area_def);
280 std::vector<fastjet::PseudoJet> inclusive_jets = sorted_by_pt(clust_pfo.inclusive_jets(5000.));
281
282 for (size_t i = 0; i < inclusive_jets.size(); i++) {
283 xAOD::Jet* jet= new xAOD::Jet();
284 xAOD::JetFourMom_t tempjetp4(inclusive_jets[i].pt(),
285 inclusive_jets[i].eta(),
286 inclusive_jets[i].phi(),
287 inclusive_jets[i].m());
288 xAOD::JetFourMom_t newArea(inclusive_jets[i].area_4vector().perp(),
289 inclusive_jets[i].area_4vector().eta(),
290 inclusive_jets[i].area_4vector().phi(),
291 inclusive_jets[i].area_4vector().m());
292 vertjets->push_back(jet);
293 jet->setJetP4(tempjetp4);
294 jet->setJetP4("JetConstitScaleMomentum",tempjetp4);
295 jet->setJetP4("JetPileupScaleMomentum",tempjetp4);
296 jet->setAttribute("ActiveArea4vec",newArea);
297 jet->setAttribute("DetectorEta",jet->eta());
298 std::vector<fastjet::PseudoJet> constituents = inclusive_jets[i].constituents();
299 float chargedpart = 0;
300 for (size_t j = 0; j < constituents.size(); j++) {
301 if (charged_pfo.count(constituents[j].user_index())>=1) {
302 chargedpart += constituents[j].perp();
303 }
304 }
305 xAOD::JetFourMom_t chargejetp4(chargedpart,inclusive_jets[i].eta(),inclusive_jets[i].phi(),inclusive_jets[i].m());
306 jet->setJetP4(m_jetchargedp4,chargejetp4);
307 }
308
309 if((m_pfoJES->modify(*vertjets)).isFailure()){
310 ATH_MSG_ERROR(" Failed to calibrate PU jet container ");
311 return pu_jets;
312 }
313
314 pu_jets.jetCont = vertjets;
315 pu_jets.jetAuxCont = vertjetsAux;
316 return pu_jets;
317 }
318
319 fastjet::PseudoJet JetForwardPFlowJvtTool::pfoToPseudoJet(const xAOD::PFO* pfo, const CP::PFO_JetMETConfig_charge& theCharge, const xAOD::Vertex *vx) const {
320 TLorentzVector pfo_p4;
321 if (CP::charged == theCharge){
322 float pweight = m_weight;
323 if( (m_wpfotool->fillWeight(*pfo,pweight)).isSuccess() ){
324 // Create a PSeudojet with the momentum of the selected IParticle
325 pfo_p4= TLorentzVector(pfo->p4().Px()*pweight,pfo->p4().Py()*pweight,pfo->p4().Pz()*pweight,pfo->e()*pweight);
326 }
327 } else if (CP::neutral == theCharge){
328 pfo_p4= pfo->GetVertexCorrectedEMFourVec(*vx);
329 }
330 fastjet::PseudoJet psj(pfo_p4);
331 // User index is used to identify the xAOD object used for the PSeudoJet
332 if (CP::charged == theCharge){
333 psj.set_user_index(pfo->index());
334 }else{
335 psj.set_user_index(-1);
336 }
337
338 return psj;
339 }
340
341 fastjet::PseudoJet JetForwardPFlowJvtTool::feToPseudoJet(const xAOD::FlowElement* fe, const CP::PFO_JetMETConfig_charge& theCharge, const xAOD::Vertex *vx) const {
342 TLorentzVector fe_p4;
343 if (CP::charged == theCharge){
344 float pweight = m_weight;
345 if( (m_wpfotool->fillWeight(*fe,pweight)).isSuccess() ){
346 // Create a Peeudojet with the momentum of the selected IParticle
347 fe_p4= TLorentzVector(fe->p4().Px()*pweight,fe->p4().Py()*pweight,fe->p4().Pz()*pweight,fe->e()*pweight);
348 }
349 } else if (CP::neutral == theCharge){
351 }
352 fastjet::PseudoJet psj(fe_p4);
353 // User index is used to identify the xAOD object used for the PseudoJet
354 if (CP::charged == theCharge){
355 psj.set_user_index(fe->index());
356 }else{
357 psj.set_user_index(-1);
358 }
359
360 return psj;
361 }
362
364 if (std::abs(jet->eta())<m_etaThresh) return false;
365 if (jet->pt()<m_forwardMinPt || (m_forwardMaxPt>0 && jet->pt()>m_forwardMaxPt) ) return false;
366 return true;
367 }
368
370 if (std::abs(jet->eta())>m_etaThresh) return false;
371 if (jet->pt()<m_centerMinPt || (m_centerMaxPt>0 && jet->pt()>m_centerMaxPt)) return false;
372 return true;
373 }
374
376 double Rpt;
377 Rpt= jet->jetP4(m_jetchargedp4).Pt()/ jet->pt();
378 return Rpt;
379 }
380
381 std::size_t JetForwardPFlowJvtTool::getPV() const{
382 if (m_includePV) return -1;
383
384 //const xAOD::VertexContainer *vxCont = 0;
386 ATH_MSG_DEBUG("Successfully retrieved primary vertex container");
387 for(const xAOD::Vertex *vx : *vxContHandle) {
388 if(vx->vertexType()==xAOD::VxType::PriVtx) return vx->index();
389 }
390 // If no verticies are found in the event the Primary Vertex container will just contain a dummy vertex and no primary vertex
391 if(vxContHandle->empty() ){
392 ATH_MSG_ERROR("Primary vertex container is empty ");
393 } else if(vxContHandle->size() != 1 ){
394 ATH_MSG_WARNING("Couldn't identify the hard-scatter primary vertex (no vertex with \"vx->vertexType()==xAOD::VxType::PriVtx\" in the container)! ");
395 }
396 // this almost certainly isn't what we should do here, the
397 // caller doesn't check this for errors
398 return 0;
399 }
400
404
405 for(const xAOD::Jet *jet : *jets) {
406 bool ishs = false;
407 bool ispu = true;
408 for(const xAOD::Jet *tjet : *truthJets) {
409 if (tjet->p4().DeltaR(jet->p4())<0.3 && tjet->pt()>10e3) ishs = true;
410 if (tjet->p4().DeltaR(jet->p4())<0.6) ispu = false;
411 }
412 isHSHandle(*jet)=ishs;
413 isPUHandle(*jet)=ispu;
414 }
415 return StatusCode::SUCCESS;
416 }
417
Scalar eta() const
pseudorapidity method
Scalar perp() const
perp method - perpendicular length
Scalar phi() const
phi method
#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.
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce(T &h)
SG::ReadDecorHandleKey< xAOD::PFO > m_orKey
Gaudi::Property< bool > m_suppressInputDependence
SG::WriteDecorHandleKey< xAOD::JetContainer > m_fjvtKey
Gaudi::Property< double > m_etaThresh
bool hasCloseByHSjet(const xAOD::Jet *jet, const xAOD::JetContainer *pjets) const
Gaudi::Property< double > m_rptCut
SG::ReadDecorHandleKey< xAOD::JetContainer > m_passJvtKey
Gaudi::Property< double > m_forwardMinPt
virtual StatusCode initialize() override
Dummy implementation of the initialisation function.
Gaudi::Property< std::string > m_jetContainerName
Gaudi::Property< double > m_centerMaxPt
fastjet::PseudoJet feToPseudoJet(const xAOD::FlowElement *fe, const CP::PFO_JetMETConfig_charge &theCharge, const xAOD::Vertex *vx) const
SG::ReadHandleKey< xAOD::VertexContainer > m_vxContKey
SG::WriteDecorHandleKey< xAOD::JetContainer > m_isHSKey
SG::ReadHandleKey< xAOD::FlowElementContainer > m_FEKey
virtual std::vector< TVector2 > calculateVertexMomenta(const xAOD::JetContainer *jets, int pvind, int vertices) const
Gaudi::Property< std::string > m_jetchargedp4
Gaudi::Property< float > m_weight
pflow::puJets buildPFlowPUjets(const xAOD::Vertex &vx) const
Gaudi::Property< double > m_centerMinPt
Gaudi::Property< bool > m_includePV
Gaudi::Property< double > m_fjvtThresh
ToolHandle< CP::WeightPFOTool > m_wpfotool
float getFJVT(const xAOD::Jet *jet, const std::vector< TVector2 > &pileupMomenta) const
SG::ReadDecorHandleKey< xAOD::FlowElement > m_orFEKey
Gaudi::Property< double > m_dzCut
ToolHandle< IJetCalibrationTool > m_pfoJES
Gaudi::Property< double > m_maxRap
SG::ReadHandleKey< xAOD::EventInfo > m_eventInfoKey
Gaudi::Property< int > m_pvind
bool isForwardJet(const xAOD::Jet *jet) const
virtual ~JetForwardPFlowJvtTool()
Destructor:
Gaudi::Property< bool > m_tightOP
bool isCentralJet(const xAOD::Jet *jet) const
StatusCode tagTruth(const xAOD::JetContainer *jets, const xAOD::JetContainer *truthJets)
SG::ReadHandleKey< xAOD::PFOContainer > m_PFOKey
Gaudi::Property< std::string > m_jetsName
SG::WriteDecorHandleKey< xAOD::JetContainer > m_isPUKey
virtual StatusCode decorate(const xAOD::JetContainer &jetCont) const override
Decorate a jet collection without otherwise modifying it.
Gaudi::Property< double > m_forwardMaxPt
Gaudi::Property< int > m_vertices
fastjet::PseudoJet pfoToPseudoJet(const xAOD::PFO *pfo, const CP::PFO_JetMETConfig_charge &theCharge, const xAOD::Vertex *vx) const
double getRpt(const xAOD::Jet *jet) const
JetForwardPFlowJvtTool(const std::string &name)
Constructor with parameters:
Gaudi::Property< double > m_neutMaxRap
SG::WriteDecorHandleKey< xAOD::JetContainer > m_fjvtRawKey
SG::ReadHandleKey< jet::TrackVertexAssociation > m_tvaKey
size_t index() const
Return the index of this element within its container.
Handle class for reading a decoration on an object.
virtual bool isValid() override final
Can the handle be successfully dereferenced?
Handle class for adding a decoration to an object.
AsgTool(const std::string &name)
Constructor specifying the tool instance's name.
Definition AsgTool.cxx:58
virtual double e() const override
The total energy of the particle.
virtual FourMom_t p4() const override
The full 4-momentum of the particle.
virtual FourMom_t p4() const
The full 4-momentum of the particle.
Definition PFO_v1.cxx:95
TLorentzVector GetVertexCorrectedEMFourVec(const xAOD::Vertex &vertexToCorrectTo) const
Correct EM scale 4-vector to point at a vertex.
Definition PFO_v1.cxx:737
virtual double e() const
The total energy of the particle.
Definition PFO_v1.cxx:81
float z() const
Returns the z position.
TLorentzVector getVertexCorrectedFourVec(const xAOD::FlowElement &fe, const xAOD::Vertex &vertexToCorrectTo)
Definition FEHelpers.cxx:13
void seedsFromEventInfo(const xAOD::EventInfo *ei, std::vector< int > &seeds)
Fill seeds vector from run & event number. This functio is separated from the class so it's easier to...
SG::ReadCondHandle< T > makeHandle(const SG::ReadCondHandleKey< T > &key, const EventContext &ctx=Gaudi::Hive::currentContext())
@ PileUp
Pile-up vertex.
@ PriVtx
Primary vertex.
Jet_v1 Jet
Definition of the current "jet version".
PFO_v1 PFO
Definition of the current "pfo version".
Definition PFO.h:17
FlowElement_v1 FlowElement
Definition of the current "pfo version".
Definition FlowElement.h:16
TrackParticle_v1 TrackParticle
Reference the current persistent version:
Vertex_v1 Vertex
Define the latest version of the vertex class.
JetContainer_v1 JetContainer
Definition of the current "jet container version".
ROOT::Math::LorentzVector< ROOT::Math::PtEtaPhiM4D< double > > JetFourMom_t
Base 4 Momentum type for Jet.
Definition JetTypes.h:17
std::shared_ptr< xAOD::JetAuxContainer > jetAuxCont
std::shared_ptr< xAOD::JetContainer > jetCont