ATLAS Offline Software
Loading...
Searching...
No Matches
PartonJetAugmentationTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
6
8#include "xAODJet/Jet.h"
11#include "fastjet/ClusterSequence.hh"
12#include "fastjet/PseudoJet.hh"
13
15#include "TLorentzVector.h"
16#include <cmath> //std::abs
17
18namespace DerivationFramework {
19 // ============================
20 // Initialize
21 // ============================
23 ATH_MSG_INFO("Initializing PartonJetAugmentationTool");
24
25 ATH_CHECK(m_truthKey.initialize());
26 ATH_CHECK(m_partonJetsKey.initialize());
27
28 return StatusCode::SUCCESS;
29 }
30
31
32 // ============================
33 // EVENT FUNCTION (runs per event)
34 // ============================
35 StatusCode PartonJetAugmentationTool::addBranches(const EventContext& ctx) const
36 {
37
38 // Get truth particles
39 //const xAOD::TruthParticleContainer* truth = nullptr;
40 //ATH_CHECK(evtStore()->retrieve(truth, "TruthParticles"));
42 if (!truth.isValid()) {
43 ATH_MSG_ERROR("Failed to retrieve TruthParticles");
44 return StatusCode::FAILURE;
45 }
46
47
48 // Create output jet container
49 auto jets = std::make_unique<xAOD::JetContainer>();
50 auto jetsAux = std::make_unique<xAOD::JetAuxContainer>();
51 jets->setStore(jetsAux.get());
52
53 TLorentzVector lj, slj, tlj;
54
55 double jet_radius = 0.4;
56 double min_pT = 20;
57
58 std::vector<TLorentzVector> ttbarDecayProducts;
59 std::vector<int> decayProduct_pdgID;
60
61 for(const auto* p : *truth){
62
63 int pdg = p->pdgId();
64
65 // Only consider W bosons and b quarks
66 if(std::abs(pdg) != 5 && std::abs(pdg) != 24)
67 continue;
68
69 bool fromTop = false;
70
71 // check parents
72 for(unsigned int ip = 0; ip < p->nParents(); ++ip){
73 const xAOD::TruthParticle* parent = p->parent(ip);
74
75 if(parent && parent->absPdgId() == 6){
76 fromTop = true;
77 break;
78 }
79 }
80
81 if(!fromTop) continue;
82 TLorentzVector vec;
83 vec.SetPtEtaPhiM(p->pt(), p->eta(), p->phi(), p->m());
84
85 ttbarDecayProducts.push_back(vec);
86 decayProduct_pdgID.push_back(pdg);
87 }
88
89 bool ok = extrajet(truth.cptr(), lj, slj, tlj, ttbarDecayProducts, decayProduct_pdgID, jet_radius, min_pT);
90
91
92 auto saveJet = [&](const TLorentzVector& v) {
93 if(v.Pt() <= 0) return; // skip empty jets
94 xAOD::Jet* jet = new xAOD::Jet();
95 jets->push_back(jet);
96 jet->setJetP4(xAOD::JetFourMom_t(v.Pt(), v.Eta(), v.Phi(), v.M()));
97 };
98
99 if(ok){
100 saveJet(lj);
101 saveJet(slj);
102 saveJet(tlj);
103 }
104
105 // Store in output DAOD
106 //ATH_CHECK(evtStore()->record(std::move(jets), "PartonJets"));
107 // ATH_CHECK(evtStore()->record(std::move(jetsAux), "PartonJetsAux."));
109 ATH_CHECK(jetsHandle.record(std::move(jets), std::move(jetsAux)));
110
111
112
113 return StatusCode::SUCCESS;
114 }
115
116
117
119
120 const xAOD::TruthParticle* parent = particle;
121 bool isLastBeforeHadronization = false;
122
123 if (particle->nChildren() > 0)
124 {
125 if (particle->child(0)->absPdgId() > 37)
126 {
128 }//if
129 if (particle->nChildren() > 1) //safety measure because some particles have only themselves as child
130 {
131 if (particle->child(1)->absPdgId() > 37)
132 {
134 }//if
135 }//if
136 }//if
137
138 // check if not after hadronization
139 while(parent->nParents()>0 && isLastBeforeHadronization == true)
140 {
141
142 parent = parent->parent(0);
143
144 if (parent->absPdgId() > 37 && parent->absPdgId() != 2212) //proton
145 {
147 }
148
149 }//while
150
152 }
153
154
156 const xAOD::TruthParticleContainer* truthParticles,
157 TLorentzVector& lj,
158 TLorentzVector& slj,
159 TLorentzVector& /*tlj unused*/,
160 const std::vector<TLorentzVector>& ttbarDecayProducts,
161 const std::vector<int>& decayProduct_pdgID,
162 double Rparam,
163 double pt_min) const {
164
165 // if true, decay products of the top quark and their children are removed from clustering
166 const bool exclude_topdecay = true;
167 std::vector<fastjet::PseudoJet> vec_status62; //partons just before top decay
168
169 for (const xAOD::TruthParticle* particle : *truthParticles) {
170 // particles with status 81-100 are for internal MC only and are not to be used in clustering
171 // particles with status > 100 are hadrons or other unstable particle
172 // They're not relevant for parton jets
173 if (particle->absPdgId() > 80 || (particle->absPdgId() > 22 && particle->absPdgId() < 38)) continue; //exclude bosons which are not gluons or photons
174
175 // Find last parton before hadronization (corresponding to status 62 in pythia)
176 //satus 11 is Pythia, status 71 is Herwig - we want to run for both
177 if((particle->status()==11 || particle->status()==71) && particle->absPdgId()!=6 ){
178
179 const xAOD::TruthParticle* topParent = particle;
180 bool exclude=false;
181
182 if (!isLastBeforeHadronization(particle))
183 {
184 continue;
185 }
186
187 //Remove all partons that come from the top quark decay
188 if (exclude_topdecay)
189 {
190 while(topParent->nParents()>0 )
191 {
192 topParent = topParent->parent(0);
193 if (topParent->nChildren() > 0 )
194 {
195 if( ( topParent->absPdgId()==24 && topParent->nParents()>0 && topParent->parent(0)->absPdgId()==6 ) ||
196 ( topParent->absPdgId()==5 && topParent->nParents()>0 && topParent->parent(0)->absPdgId()==6 ) ) exclude=true;
197 }
198 if (exclude) break;
199 }//while
200 }// if (exclude_topdecay)
201
202
203 if(!exclude){
204 fastjet::PseudoJet tmp(particle->px(), particle->py(), particle->pz(), particle->e());
205 vec_status62.push_back(tmp);
206 } //if(!exclude)
207
208 }
209 }
210
211 // Fastjet analysis - select algorithm and parameters
212 fastjet::JetDefinition jetDef(fastjet::antikt_algorithm, Rparam, fastjet::E_scheme, fastjet::Best);
213
214 std::vector <fastjet::PseudoJet> inclusivePertJet62, sortedPertJet62;
215
216 fastjet::ClusterSequence clustSeqPertJet62(vec_status62, jetDef);
217 // Extract inclusive jets sorted by pT, and ask for Pt>pt_min
218 inclusivePertJet62 = clustSeqPertJet62.inclusive_jets(pt_min * 1000.); //GeV
219 sortedPertJet62 = fastjet::sorted_by_pt(inclusivePertJet62);
220 // std::cout<<"containers size: "<<inclusivePertJet22.size()<<" "<<inclusivePertJet22.size()<<"\n";
221
222 if(inclusivePertJet62.size()==0){
223 // std::cout<<"Extrajet containers are empty! No partonic extrajet!\n";
224 return false;
225 }
226
227 if (inclusivePertJet62.size()!=0)
228 {
229 // remove jets that are in deltaR <= 0.4 to one of the top decay products
230 for (long unsigned int i = 0; i < ttbarDecayProducts.size(); ++i)
231 {
232 for (long unsigned int j = 0; j < sortedPertJet62.size(); ++j)
233 {
234 if (sortedPertJet62.at(j).delta_R(ttbarDecayProducts[i]) <= 0.4
235 && std::abs(decayProduct_pdgID[i]) != 11 // electron
236 && std::abs(decayProduct_pdgID[i]) != 12 // electron neutrino
237 && std::abs(decayProduct_pdgID[i]) != 13 // muon
238 && std::abs(decayProduct_pdgID[i]) != 14 // muon neutrino
239 // do not veto tau leptons since they also decay hadronically
240 && std::abs(decayProduct_pdgID[i]) != 16 // tau neutrino
241 )
242 {
243
244 sortedPertJet62.erase(sortedPertJet62.begin() + j);
245
246 --j;
247 }
248 }//for
249
250 }//for
251
252 bool leading_jet_found = false;
253
254
255 /* The following variables are now unused, left for review
256 bool third_jet_found = false;
257 bool subleading_jet_found = false;
258 */
259
260 for (std::vector<fastjet::PseudoJet>::iterator jet = sortedPertJet62.begin(); jet != sortedPertJet62.end(); ++jet)
261 {
262 if(!leading_jet_found){
263 lj.SetPxPyPzE(jet->px(),jet->py(),jet->pz(),jet->e()); //only saves leading jet
264 leading_jet_found = true;
265 }
266 else { //leading_jet_found is 'true' and subleading_jet_found is 'false'
267 slj.SetPxPyPzE(jet->px(),jet->py(),jet->pz(),jet->e()); //only saves leading jet
268 /* the following is set but then never used
269 subleading_jet_found = true;
270 */
271 break;//exits loop if leading_jet_found
272 }
280 } //for
281 }
282
283
284 if (lj.Pt() >= pt_min*1000.) return true;
285
286 return false;
287
288 }
289
290}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
std::vector< size_t > vec
Helpers for checking error return status codes and reporting errors.
SG::WriteHandleKey< xAOD::JetContainer > m_partonJetsKey
SG::ReadHandleKey< xAOD::TruthParticleContainer > m_truthKey
bool isLastBeforeHadronization(const xAOD::TruthParticle *p) const
bool extrajet(const xAOD::TruthParticleContainer *truthParticles, TLorentzVector &lj, TLorentzVector &slj, TLorentzVector &tlj, const std::vector< TLorentzVector > &ttbarDecayProducts, const std::vector< int > &decayProduct_pdgID, double Rparam, double pt_min) const
virtual StatusCode addBranches(const EventContext &ctx) const override
virtual bool isValid() override final
Can the handle be successfully dereferenced?
const_pointer_type cptr()
Dereference the pointer.
StatusCode record(std::unique_ptr< T > data)
Record a const object to the store.
const TruthParticle_v1 * parent(size_t i) const
Retrieve the i-th mother (TruthParticle) of this TruthParticle.
int absPdgId() const
Absolute PDG ID code (often useful).
size_t nParents() const
Number of parents of this particle.
size_t nChildren() const
Number of children of this particle.
std::set< std::string > exclude
list of directories to be excluded
Definition hcg.cxx:100
THE reconstruction tool.
Jet_v1 Jet
Definition of the current "jet version".
TruthParticle_v1 TruthParticle
Typedef to implementation.
TruthParticleContainer_v1 TruthParticleContainer
Declare the latest version of the truth particle container.
ROOT::Math::LorentzVector< ROOT::Math::PtEtaPhiM4D< double > > JetFourMom_t
Base 4 Momentum type for Jet.
Definition JetTypes.h:17