ATLAS Offline Software
HadronOriginClassifier.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3 */
6 #include "StoreGate/ReadHandle.h"
7 
8 namespace {
10  struct Sample {
12 
14  Sample(int low, int high, GEN_id gen) :
15  low(low), high(high), gen(gen) {}
16 
18  Sample(int id, GEN_id gen) :
19  Sample(id, id, gen) {}
20 
21  int low{};
22  int high{};
23  GEN_id gen{GEN_id::Pythia6};
24  };
25 }
26 
27 namespace DerivationFramework{
28 
29  HadronOriginClassifier::HadronOriginClassifier(const std::string& t, const std::string& n, const IInterface* p):
30  AthAlgTool(t,n,p)
31  {
32  }
33 
35 
37  ATH_MSG_INFO("Initialize " );
38  ATH_MSG_INFO("DSID " << m_DSID );
39 
40  ATH_CHECK(m_mcName.initialize());
41 
42  static const std::vector<Sample> samples = {
43  // all Herwig++/Herwig7 showered samples
44  {346346, 346348, GEN_id::HerwigPP},
45  {410003, GEN_id::HerwigPP}, {410008, GEN_id::HerwigPP}, //aMC@NLO+Hpp
46  {410004, GEN_id::HerwigPP}, {410163, GEN_id::HerwigPP}, //Powheg+Hpp
47  {410232, 410233, GEN_id::HerwigPP}, //first attempt for Powheg+H7 / aMC@NLO+H7
48  {410525, 410530, GEN_id::HerwigPP}, //New Powheg+H7 samples
49  {407037, 407040, GEN_id::HerwigPP}, //Powheg+Hpp MET/HT sliced
50  {410536, 410537, GEN_id::HerwigPP}, {410245, GEN_id::HerwigPP}, //aMC@NLO+H++ , ttbb
51  {410557, 410559, GEN_id::HerwigPP}, // new Powheg+H7, mc16
52  {411082, 411090, GEN_id::HerwigPP}, //Powheg+H7 HF-filtered
53  {407354, 407356, GEN_id::HerwigPP}, //Powheg+H7 ttbar HT-filtered
54  {411233, 411234, GEN_id::HerwigPP}, //Powheg+H7.1.3 ttbar
55  {411316, GEN_id::HerwigPP}, //Powheg+H7 allhad ttbar
56  {411329, 411334, GEN_id::HerwigPP}, //Powheg+H7.1.3 ttbar HF-filtered
57  {411335, 411337, GEN_id::HerwigPP}, //Powheg+H7.1.3 ttbar HT-filtered
58  {412116, 412117, GEN_id::HerwigPP}, //amc@NLO+H7.1.3 ttbar
59  {504329, GEN_id::HerwigPP}, {504333, GEN_id::HerwigPP}, {504341, GEN_id::HerwigPP}, //amc@NLO+H7.2.1 refined ttZ
60  {601239, 601240, GEN_id::HerwigPP},
61  {601668, GEN_id::HerwigPP},
62  {603905, 603906, GEN_id::HerwigPP}, // ttbb Powheg+H7 dilep, ljet, allhad
63  {504337, GEN_id::HerwigPP}, {504345, GEN_id::HerwigPP}, // aMC@NLO+H7 ttZ
64  {526034, GEN_id::HerwigPP}, // aMC@NLO+H7 4tops
65  {600666, 600667, GEN_id::HerwigPP}, // Powheg+H7 ttbar H7UE
66  {601414, 601415, GEN_id::HerwigPP}, // Powheg+H7 ttbar A14
67  {602635, 602635, GEN_id::HerwigPP}, // Powheg+H7 ttH PDF4LHC21
68  {602846, 602849, GEN_id::HerwigPP}, // Powheg+H7 ttW NNPDF30NLO EW
69 
70  // all Pythia8 showered samples
71  {304014, GEN_id::Pythia8}, // amc@NLO+P8 3top
72  {346229, 346234, GEN_id::Pythia8}, // amc@NLO+P8 tHjb
73  {410006, GEN_id::Pythia8}, //Powheg+P8 old main31
74  {410081, GEN_id::Pythia8}, //amc@NLO+P8 ttV
75  {410500, GEN_id::Pythia8}, //Powheg+P8 new main31, hdamp=mt
76  {410501, 410508, GEN_id::Pythia8}, //Powheg+P8 new main31, hdamp=1.5m // Boosted samples are included 410507 410508
77  {410511, 410524, GEN_id::Pythia8}, //Powheg+P8 new main31, hdamp=1.5mt, radiation systematics
78  {410531, 410535, GEN_id::Pythia8}, //Powheg+P8 allhad samples
79  {346343, 346345, GEN_id::Pythia8}, //Powheg+P8 ttH
80  {412123, GEN_id::Pythia8}, // MG+P8 ttW
81  {410155, GEN_id::Pythia8}, // aMC@NlO+P8 ttW
82  {410159, 410160, GEN_id::Pythia8}, //aMC@NLO+P8, old settings
83  {410218, 410220, GEN_id::Pythia8}, // aMC@NlO+P8 ttZ
84  {410276, 410278, GEN_id::Pythia8}, // aMC@NlO+P8 ttZ_lowMass
85  {410225, 410227, GEN_id::Pythia8}, {410274, 410275, GEN_id::Pythia8}, //aMC@NLO+P8, new settings
86  {410568, 410569, GEN_id::Pythia8}, // nonallhad boosted c-filtered
87  {410244, GEN_id::Pythia8}, //aMC@NLO+P8, ttbb (old)
88  {410441, 410442, GEN_id::Pythia8}, //new aMC@NLO+P8 mc16, new shower starting scale
89  {410464, 410466, GEN_id::Pythia8}, //new aMC@NLO+P8 mc16, new shower starting scale, no shower weights
90  {410470, 410472, GEN_id::Pythia8}, {410480, 410482, GEN_id::Pythia8}, //new Powheg+P8 mc16
91  {410452, GEN_id::Pythia8}, //new aMC@NLO+P8 FxFx mc16
92  {411073, 411081, GEN_id::Pythia8}, //Powheg+P8 HF-filtered
93  {412043, 412044, GEN_id::Pythia8}, {500326, GEN_id::Pythia8}, //aMC@NLO+P8 4top
94  {412066, 412074, GEN_id::Pythia8}, //aMC@NLO+P8 HF-filtered
95  {411068, 411070, GEN_id::Pythia8}, //Powheg+P8 ttbb
96  {410265, 410267, GEN_id::Pythia8}, //aMC@NLO+P8 ttbb
97  {411178, 411180, GEN_id::Pythia8}, {411275, GEN_id::Pythia8}, //Powheg+P8 ttbb OTF production - ATLMCPROD-7240
98  {501720, GEN_id::Pythia8}, // aMC@NLO+P8 FxFx ttW
99  {600791, 600792, GEN_id::Pythia8}, //Powheg+P8 ttbb - ATLMCPROD-9179
100  {600737, 600738, GEN_id::Pythia8}, //Powheg+P8 ttbb - ATLMCPROD-9179
101  {601226, 601228, GEN_id::Pythia8}, // Powheg+P8 ttbb bornzerodamp cut 5, ATLMCPROD-9694
102  {407342, 407344, GEN_id::Pythia8}, {411391, GEN_id::Pythia8}, //Powheg+P8 ttbar HT-filtered
103  {407345, 407347, GEN_id::Pythia8}, //Powheg+P8 ttbar MET-filtered
104  {407348, 407350, GEN_id::Pythia8}, //aMC@NLO+P8 ttbar HT-filtered
105  {504330, 504332, GEN_id::Pythia8}, {504334, 504336, GEN_id::Pythia8}, {504338, GEN_id::Pythia8}, {504342, 504344, GEN_id::Pythia8}, {504346, GEN_id::Pythia8}, //aMC@NLO+P8 refined ttZ
106  {601491, 601492, GEN_id::Pythia8}, //Pow+Py8 ttbar pTHard variations - ATLMCPROD-10168
107  {601495, 601498, GEN_id::Pythia8}, //Pow+Py8 ttbar pTHard variations - ATLMCPROD-10168
108  {601229, 601230, GEN_id::Pythia8}, // mc23 ttbar dilep, singlelep
109  {601237, GEN_id::Pythia8}, // mc23 ttbar allhad
110  {601398, 601399, GEN_id::Pythia8}, // mc23 ttbar dilep, singlelep hdamp517p5
111  {601491, GEN_id::Pythia8}, {601495, GEN_id::Pythia8}, {601497, GEN_id::Pythia8}, // mc23 ttbar pThard variations, dilep, singlelep, allhad
112  {601783, 601784, GEN_id::Pythia8}, // Powheg+P8 ttbb bornzerodamp cut 5 pThard variations - ATLMCPROD-10527
113  {603003, 603004, GEN_id::Pythia8}, // Powheg+P8 ttbb nominal and pthard1 allhad
114  {603190, 603193, GEN_id::Pythia8}, // Powheg+P8 ttbb nominal and pthard1 dilep, ljet
115  {604482, 604483, GEN_id::Pythia8}, // mc23 Powheg+P8 ttbar recoilToTop
116  {411369, 411374, GEN_id::Pythia8}, // PP8 ttbar Var1
117  {411290, GEN_id::Pythia8}, // PP8 ttbar rb1p05
118  {508792, GEN_id::Pythia8}, // aMC@NLO+P8 ttbar smeftsim
119  {510203, GEN_id::Pythia8}, // aMC@NLO+P8 4tops
120  {522035, 522038, GEN_id::Pythia8}, // aMC@NLO+P8 ttZqq
121  {523243, GEN_id::Pythia8}, // aMC@NLO+P8 4tops
122  {601356, 601357, GEN_id::Pythia8}, // PP8 ttbar Trec
123  {602067, 602072, GEN_id::Pythia8}, // PP8 ttH pthard1/2
124  {602637, GEN_id::Pythia8}, // PP8 ttH PDF4LHC21
125  {602852, 602853, GEN_id::Pythia8}, // PP8 ttH PDF4LHC21 pthard1
126  {602886, 602889, GEN_id::Pythia8}, // PP8 ttW NNPDF23
127  {603851, 603854, GEN_id::Pythia8}, // PP8 ggH
128  {604224, GEN_id::Pythia8}, // PP8 ttH allhad HTop
129 
130  // all Sherpa showered samples
131  {410186, 410189, GEN_id::Sherpa}, //Sherpa 2.2.0
132  {410249, 410252, GEN_id::Sherpa}, //Sherpa 2.2.1
133  {410342, 410347, GEN_id::Sherpa}, //Sherpa 2.2.1 sys
134  {410350, 410355, GEN_id::Sherpa}, //Sherpa 2.2.1 sys
135  {410357, 410359, GEN_id::Sherpa}, //Sherpa 2.2.1 sys
136  {410361, 410367, GEN_id::Sherpa}, //Sherpa 2.2.1 sys
137  {410281, 410283, GEN_id::Sherpa}, //Sherpa BFilter
138  {410051, GEN_id::Sherpa}, //Sherpa ttbb (ICHEP sample)
139  {410323, 410325, GEN_id::Sherpa}, {410369, GEN_id::Sherpa}, //New Sherpa 2.2.1 ttbb
140  {364345, 364348, GEN_id::Sherpa}, //Sherpa 2.2.4 (test)
141  {410424, 410427, GEN_id::Sherpa}, //Sherpa 2.2.4
142  {410661, 410664, GEN_id::Sherpa}, //Sherpa 2.2.4 ttbb
143  {421152, 421158, GEN_id::Sherpa}, //Sherpa2.2.8 ttbar
144  {413023, GEN_id::Sherpa}, // sherpa 2.2.1 ttZ
145  {700000, GEN_id::Sherpa}, // Sherpa 2.2.8 ttW
146  {700168, GEN_id::Sherpa}, // Sherpa 2.2.10 ttW
147  {700205, GEN_id::Sherpa}, // Sherpa 2.2.10 ttW EWK
148  {700309, GEN_id::Sherpa}, // Sherpa 2.2.11 ttZ
149  {700051, 700054, GEN_id::Sherpa}, //Sherpa2.2.8 ttbb
150  {700121, 700124, GEN_id::Sherpa}, //Sherpa2.2.10 ttbar
151  {700164, 700167, GEN_id::Sherpa}, //Sherpa2.2.10 ttbb
152  {700807, 700809, GEN_id::Sherpa}, //Sherpa2.2.14 ttbar
153  {700659, 700662, GEN_id::Sherpa}, // Sherpa 2.2.12 ttbar maxHTavrgTopPT
154  {700712, GEN_id::Sherpa}, // Sherpa 2.2.14 4tops muQHT2
155  {700986, 700997, GEN_id::Sherpa}, // Sherpa 2.2.14 ttW
156  {701251, 701259, GEN_id::Sherpa}, // Sherpa ttW
157 
158  // everything else from HFDSIDList
159  // and the mc20/mc23 central pages as of Sept 2025
160  {301528, 301532, GEN_id::Pythia6},
161  {301539, GEN_id::Sherpa},
162  {302910, 302924, GEN_id::Sherpa},
163  {303480, 303487, GEN_id::Sherpa},
164  {303722, 303726, GEN_id::Pythia6},
165  {306600, 306617, GEN_id::Sherpa},
166  {307479, 307502, GEN_id::Sherpa},
167  {343362, GEN_id::Pythia6},
168  {343431, 343434, GEN_id::Pythia8},
169  {343637, GEN_id::Pythia6},
170  {343852, 343854, GEN_id::Pythia8},
171  {344171, GEN_id::Pythia6},
172  {345935, GEN_id::Pythia8},
173  {345951, GEN_id::Pythia8},
174  {346031, GEN_id::Pythia8},
175  {407009, 407012, GEN_id::Pythia6},
176  {407029, 407036, GEN_id::Pythia6},
177  {407041, 407048, GEN_id::Pythia8},
178  {407200, 407204, GEN_id::Pythia8},
179  {407320, 407321, GEN_id::Pythia8},
180  {407322, 407323, GEN_id::Pythia6},
181  {407324, 407335, GEN_id::Pythia8},
182  {407351, 407353, GEN_id::Pythia8},
183  {407357, 407359, GEN_id::HerwigPP},
184  {410000, 410002, GEN_id::Pythia6},
185  {410007, GEN_id::Pythia6},
186  {410009, GEN_id::Pythia6},
187  {410021, 410024, GEN_id::Sherpa},
188  {410028, 410029, GEN_id::Pythia8},
189  {410037, 410046, GEN_id::Pythia6},
190  {410052, 410058, GEN_id::Sherpa},
191  {410060, GEN_id::HerwigPP},
192  {410066, 410070, GEN_id::Pythia8},
193  {410073, GEN_id::Pythia8},
194  {410075, GEN_id::Pythia8},
195  {410082, 410084, GEN_id::Pythia8},
196  {410087, 410089, GEN_id::Pythia8},
197  {410111, 410116, GEN_id::Pythia8},
198  {410120, 410121, GEN_id::Pythia6},
199  {410142, 410144, GEN_id::Sherpa},
200  {410156, 410157, GEN_id::Pythia8},
201  {410161, 410162, GEN_id::Pythia6},
202  {410175, 410176, GEN_id::Pythia8},
203  {410178, 410179, GEN_id::Pythia8},
204  {410181, 410182, GEN_id::Pythia8},
205  {410184, 410185, GEN_id::Sherpa},
206  {410190, 410213, GEN_id::Pythia8},
207  {410234, GEN_id::Pythia8},
208  {410248, GEN_id::Pythia8},
209  {410257, 410262, GEN_id::Pythia8},
210  {410284, 410322, GEN_id::Pythia8},
211  {410326, 410339, GEN_id::Sherpa},
212  {410368, GEN_id::Pythia8},
213  {410370, 410394, GEN_id::Pythia8},
214  {410395, 410396, GEN_id::HerwigPP},
215  {410397, 410399, GEN_id::Pythia8},
216  {410404, 410405, GEN_id::Pythia8},
217  {410410, 410413, GEN_id::Pythia8},
218  {410420, 410423, GEN_id::Sherpa},
219  {410428, 410431, GEN_id::Pythia8},
220  {410432, 410433, GEN_id::HerwigPP},
221  {410434, 410435, GEN_id::Pythia8},
222  {410444, 410445, GEN_id::Pythia8},
223  {410446, GEN_id::HerwigPP},
224  {410447, GEN_id::Pythia8},
225  {410467, GEN_id::Pythia8},
226  {410468, GEN_id::HerwigPP},
227  {410469, GEN_id::Pythia8},
228  {410491, 410498, GEN_id::Pythia8},
229  {410509, GEN_id::Pythia8},
230  {410542, GEN_id::Pythia8},
231  {410544, GEN_id::Pythia8},
232  {410545, 410546, GEN_id::HerwigPP},
233  {410555, GEN_id::HerwigPP},
234  {410633, 410637, GEN_id::Pythia8},
235  {410690, 410691, GEN_id::Pythia8},
236  {411000, 411005, GEN_id::Pythia8},
237  {411044, 411059, GEN_id::Pythia8},
238  {411125, 411142, GEN_id::HerwigPP},
239  {411143, 411160, GEN_id::Pythia8},
240  {411168, GEN_id::HerwigPP},
241  {411169, GEN_id::Pythia8},
242  {411235, 411268, GEN_id::Pythia8},
243  {411269, GEN_id::HerwigPP},
244  {411278, GEN_id::Pythia6},
245  {411279, 411281, GEN_id::Pythia8},
246  {411282, GEN_id::HerwigPP},
247  {411286, GEN_id::Pythia8},
248  {411288, 411289, GEN_id::Pythia8},
249  {411292, 411307, GEN_id::Pythia8},
250  {411310, GEN_id::Pythia8},
251  {412008, 412014, GEN_id::Pythia8},
252  {412017, 412025, GEN_id::Pythia8},
253  {412028, 412036, GEN_id::Pythia8},
254  {412039, 412040, GEN_id::Pythia8},
255  {412090, 412092, GEN_id::HerwigPP},
256  {412112, 412114, GEN_id::Pythia8},
257  {412175, GEN_id::HerwigPP},
258  {413005, GEN_id::Sherpa},
259  {413008, GEN_id::Sherpa},
260  {413022, GEN_id::Sherpa},
261  {426072, GEN_id::Pythia6},
262  {426075, 426078, GEN_id::Pythia6},
263  {426082, GEN_id::Pythia6},
264  {426085, 426088, GEN_id::Pythia6},
265  {426090, 426097, GEN_id::Pythia6},
266  {500462, 500463, GEN_id::Pythia8},
267  {500800, GEN_id::Pythia8},
268  {502957, 502958, GEN_id::Pythia8},
269  {504553, GEN_id::HerwigPP},
270  {504554, GEN_id::Pythia8},
271  {504689, GEN_id::HerwigPP},
272  {508772, 508773, GEN_id::Pythia8},
273  {508780, GEN_id::HerwigPP},
274  {508781, GEN_id::Pythia8},
275  {508985, 508986, GEN_id::Pythia8},
276  {510212, 510213, GEN_id::Pythia8},
277  {521379, 521380, GEN_id::Pythia8},
278  {521384, 521385, GEN_id::Pythia8},
279  {522023, 522026, GEN_id::Pythia8},
280  {522027, GEN_id::HerwigPP},
281  {522028, 522030, GEN_id::Pythia8},
282  {522031, GEN_id::HerwigPP},
283  {522032, 522034, GEN_id::Pythia8},
284  {522039, GEN_id::HerwigPP},
285  {522040, 522042, GEN_id::Pythia8},
286  {542859, 542860, GEN_id::Pythia8},
287  {542867, 542868, GEN_id::Pythia8},
288  {545023, 545024, GEN_id::Pythia8},
289  {545025, 545026, GEN_id::HerwigPP},
290  {545790, 545791, GEN_id::Pythia8},
291  {561980, GEN_id::Pythia8},
292  {561982, 561990, GEN_id::Pythia8},
293  {600031, GEN_id::Pythia8},
294  {600638, 600639, GEN_id::Pythia8},
295  {600668, GEN_id::HerwigPP},
296  {600787, 600790, GEN_id::HerwigPP},
297  {600793, 600796, GEN_id::Pythia8},
298  {601284, GEN_id::Pythia8},
299  {601403, GEN_id::Pythia8},
300  {601407, GEN_id::Pythia8},
301  {601607, GEN_id::HerwigPP},
302  {601656, GEN_id::Pythia8},
303  {601669, 601670, GEN_id::Pythia8},
304  {601672, GEN_id::Pythia8},
305  {601708, GEN_id::Pythia8},
306  {602423, GEN_id::Pythia8},
307  {602636, GEN_id::HerwigPP},
308  {602638, GEN_id::Pythia8},
309  {602646, 602647, GEN_id::Pythia8},
310  {602687, 602688, GEN_id::Pythia8},
311  {602843, 602844, GEN_id::Pythia8},
312  {602850, 602851, GEN_id::Pythia8},
313  {603011, 603013, GEN_id::Pythia8},
314  {603855, 603856, GEN_id::HerwigPP},
315  {603872, 603873, GEN_id::Pythia8},
316  {603986, 603989, GEN_id::Pythia8},
317  {604008, 604009, GEN_id::HerwigPP},
318  {604018, 604021, GEN_id::Pythia8},
319  {604468, 604481, GEN_id::Pythia8},
320  {700706, GEN_id::Sherpa},
321  {700737, 700756, GEN_id::Sherpa},
322  {700810, GEN_id::Sherpa},
323  {701260, 701262, GEN_id::Sherpa},
324  {701265, 701282, GEN_id::Sherpa},
325  {802380, 802381, GEN_id::Pythia8},
326  {950787, GEN_id::Pythia8},
327  };
328 
329  // Linear search for sample and assign properties:
330  for (const auto& s : samples) {
331  if (m_DSID>=s.low && m_DSID<=s.high) {
332  m_GenUsed = s.gen;
333  return StatusCode::SUCCESS;
334  }
335  }
336 
337  // the default is Pythia6, so no need to list the Pythia6 showered samples
338  // these are:
339  // 410000-410002
340  // 410007, 410009, 410120-410121
341  // 301528-301532
342  // 303722-303726
343  // 407009-407012
344  // 407029-407036
345  // 410120
346  // 426090-426097
347  // 429007
349 
350  return StatusCode::SUCCESS;
351  }
352 
353  /*
354  --------------------------------------------------------------------------------------------------------------------------------------
355  ------------------------------------------------------------- Hadron Map -------------------------------------------------------------
356  --------------------------------------------------------------------------------------------------------------------------------------
357  */
358 
359  // Define the function GetOriginMap that determines the origin of the hadrons.
360  std::map<const xAOD::TruthParticle*, DerivationFramework::HadronOriginClassifier::HF_id> HadronOriginClassifier::GetOriginMap() const {
361  // Create a set of maps to store the information about the hadrons and the partons
362  std::map<const xAOD::TruthParticle*, int> mainHadronMap; // Map with main hadrons and their flavor.
363  std::map<const xAOD::TruthParticle*, HF_id> partonsOrigin; // Map with partons and their category (from top, W, H, MPI, FSR, extra).
364  std::map<const xAOD::TruthParticle*, const xAOD::TruthParticle*> hadronsPartons; // Map with hadrons and their matched parton.
365  std::map<const xAOD::TruthParticle*, HF_id> hadronsOrigin; // Map with hadrons and their category (from top, W, H, MPI, FSR, extra)
366  // Fill the maps mainHadronMap and partonsOrigin
367  buildPartonsHadronsMaps(mainHadronMap, partonsOrigin);
368  // Create two maps to know which partons and hadrons have already been matched.
369  std::vector<const xAOD::TruthParticle*> matched_partons;
370  std::vector<const xAOD::TruthParticle*> matched_hadrons;
371  // Use a while to go through the HF hadrons in mainHadronMap and partons in partonsOrigin.
372  while (matched_partons.size()<partonsOrigin.size() && matched_hadrons.size()<mainHadronMap.size()){
373  // Create a float variable to store the DeltaR between a parton and the closest hadron.
374  float dR=999.;
375  // Create two pointers for TruthParticle type to go through the partons and hadrons.
376  const xAOD::TruthParticle* hadron=nullptr;
377  const xAOD::TruthParticle* parton=nullptr;
378  // Use a for to go through the partonsOrigin.
379  for(std::map<const xAOD::TruthParticle*, HF_id>::iterator itr = partonsOrigin.begin(); itr!=partonsOrigin.end(); ++itr){
380  // Check if the parton has already been matched to an hadron.
381  if(std::find(matched_partons.begin(), matched_partons.end(), (*itr).first) != matched_partons.end()) continue;
382  // Extract the pt of the parton.
383  TVector3 v, vtmp;
384  if ((*itr).first->pt()>0.)
385  v.SetPtEtaPhi((*itr).first->pt(),(*itr).first->eta(),(*itr).first->phi());
386  else // Protection against FPE from eta and phi calculation
387  v.SetXYZ(0.,0.,(*itr).first->pz());
388  // Use a for to go through the HF hadrons in mainHadronMap.
389  for(std::map<const xAOD::TruthParticle*, int>::iterator it = mainHadronMap.begin(); it!=mainHadronMap.end(); ++it){
390  // Check if the hadron has already been matched to a parton.
391  if(std::find(matched_hadrons.begin(), matched_hadrons.end(), (*it).first) != matched_hadrons.end()) continue;
392  // Check if the hadron's flavour matches the one of the parton.
393  if((*it).second != (*itr).first->absPdgId()) continue;
394  // Extract the pt of the hadron.
395  vtmp.SetPtEtaPhi((*it).first->pt(),(*it).first->eta(),(*it).first->phi());
396  // Compute Delta R between hadron and parton and store in dR if it is smaller than the current value.
397  // Also store the parton and hadron in the pointers that have been previous created.
398  if(vtmp.DeltaR(v) < dR){
399  dR = vtmp.DeltaR(v);
400  hadron = (*it).first;
401  parton = (*itr).first;
402  }
403  }//loop hadrons
404  }//loop partons
405  // Add the matched part-hadron pair in the corresponding maps.
406  matched_partons.push_back(parton);
407  matched_hadrons.push_back(hadron);
408  hadronsPartons[ hadron ] = parton;
409  }
410 
411  // Use a for to go through the HF hadrons in mainHadronMap.
412  for(std::map<const xAOD::TruthParticle*, int>::iterator it = mainHadronMap.begin(); it!=mainHadronMap.end(); ++it){
413  // Extract the current hadron.
414  const xAOD::TruthParticle* hadron = (*it).first;
415  // Check if the hadron has been matched to a parton.
416  // If it has been matched to any hadron, use it to determine the origin.
417  // Otherwise, the hadron is considered extra.
418  if(hadronsPartons.find(hadron)!=hadronsPartons.end()){
419  hadronsOrigin[hadron] = partonsOrigin[ hadronsPartons[hadron] ];
420  } else{
421  hadronsOrigin[hadron] = extrajet;
422  }
423  }
424  return hadronsOrigin;
425  }
426 
427  // Define the function buildPartonsHadronsMaps that determines the flavour of the hadrons and the origin of the partons.
428  void HadronOriginClassifier::buildPartonsHadronsMaps(std::map<const xAOD::TruthParticle*,int>& mainHadronMap, std::map<const xAOD::TruthParticle*,HF_id>& partonsOrigin) const {
429  // Extract the TruthParticles container.
430  const EventContext& ctx = Gaudi::Hive::currentContext();
431  SG::ReadHandle<xAOD::TruthEventContainer> xTruthEventContainer(m_mcName, ctx);
432  if (!xTruthEventContainer.isValid()) {
433  ATH_MSG_WARNING("Could not retrieve " <<m_mcName);
434  }
435 
436  // Create a container with TruthParticles to store the hadrons that has already been saved.
437  std::set<const xAOD::TruthParticle*> usedHadron;
438  for ( const auto* truthevent : *xTruthEventContainer ) {
439  // Use a for to go through the TruthParticles.
440  for(unsigned int i = 0; i < truthevent->nTruthParticles(); i++){
441  // Extract the i-th particle.
442  const xAOD::TruthParticle* part = truthevent->truthParticle(i);
443  if(!part) continue;
444  // Simulated particles are not considered.
446  // Create a set of boolean variables to indicate the type of particle.
447  bool isbquark = false; // The particle is a b-quark.
448  bool iscquark = false; // The particle is a c-quark.
449  bool isHFhadron = false; // The particle is a HF hadron.
450  // Extract the pdgid of the particle and use it to determine the type of particle.
451  int pdgid = part->absPdgId();
452  if( MC::isBottom(pdgid) ){
453  isbquark=true;
454  }
455  else if( MC::isCharm(pdgid) ){
456  iscquark=true;
457  }
459  isHFhadron=true;
460  }
461  else{
462  continue;
463  }
464  // For HF quarks (b or c), check their category.
465  // The category is determined looking for the parents.
466  if(isbquark){
467  // In this case, the parton is a b-quark.
468  // Check the category of the b-quark.
470  partonsOrigin[ part ] = b_from_W;
471  }
472  else if(isDirectlyFromTop(part)){
473  partonsOrigin[ part ] = b_from_top;
474  }
475  else if((IsHerwigPP()||IsSherpa())&&isDirectlyFSR(part)){
476  partonsOrigin[ part ] = b_FSR;
477  }
478  else if(IsPythia8()&&isDirectlyFSRPythia8(part)){
479  partonsOrigin[ part ] = b_FSR;
480  }
481  else if(IsPythia6()&&isDirectlyFSRPythia6(part)){
482  partonsOrigin[ part ] = b_FSR;
483  }
484  else if(IsPythia6()&&isDirectlyMPIPythia6(part)){
485  partonsOrigin[ part ] = b_MPI;
486  }
487  else if(IsPythia8()&&isDirectlyMPIPythia8(part)){
488  partonsOrigin[ part ] = b_MPI;
489  }
490  else if(IsSherpa()&&isDirectlyMPISherpa(part)){
491  partonsOrigin[ part ] = b_MPI;
492  }
493  }
494  if(iscquark){
495  // In this case, the parton is a c-quark.
496  // Check the category of the b-quark.
498  partonsOrigin[ part ] = c_from_W;
499  }
500  else if(isDirectlyFromTop(part)){
501  partonsOrigin[ part ] = c_from_top;
502  }
503  else if((IsHerwigPP()&&IsSherpa())&&isDirectlyFSR(part)){
504  partonsOrigin[ part ] = c_FSR;
505  }
506  else if(IsPythia8()&&isDirectlyFSRPythia8(part)){
507  partonsOrigin[ part ] = c_FSR;
508  }
509  else if(IsPythia6()&&isDirectlyFSRPythia6(part)){
510  partonsOrigin[ part ] = c_FSR;
511  }
512  else if(IsPythia6()&&isDirectlyMPIPythia6(part)){
513  partonsOrigin[ part ] = c_MPI;
514  }
515  else if(IsPythia8()&&isDirectlyMPIPythia8(part)){
516  partonsOrigin[ part ] = c_MPI;
517  }
518  else if(IsSherpa()&&isDirectlyMPISherpa(part)){
519  partonsOrigin[ part ] = c_MPI;
520  }
521  }
522  // The HF hadrons are stored in the map mainHadronMap if they are not repeated.
523  if(isHFhadron && !isCHadronFromB(part)){
524  // In this case, the particle is a HF hadron but not a C-Hadron from a B-hadron.
525  // If the hadron is not in usedHadron, then add it in mainHadronMap with fillHadronMap function.
526  if(usedHadron.insert(part).second) {
527  fillHadronMap(usedHadron, mainHadronMap,part,part);
528  }
529  }
530  }//loop on particles
531  }//loop on truthevent container
532  }
533 
534  /*
535  ---------------------------------------------------------------------------------------------------------------------------------------
536  ------------------------------------------------------------ Particle Type ------------------------------------------------------------
537  ---------------------------------------------------------------------------------------------------------------------------------------
538  */
539  bool HadronOriginClassifier::isCHadronFromB(const xAOD::TruthParticle* part, std::shared_ptr<std::set<const xAOD::TruthParticle*>> checked ) const{
540  if(!MC::isCharmHadron(part)) return false;
541  if (!checked) checked = std::make_shared<std::set<const xAOD::TruthParticle*>>();
542  checked ->insert(part);
543 
544  for(unsigned int i=0; i<part->nParents(); ++i){
545  const xAOD::TruthParticle* parent = part->parent(i);
546  if(!parent) continue;
547  if(checked->count(parent)) continue;
548  checked->insert(parent);
549  if( MC::isBottomHadron(parent) ){
550  return true;
551  }
553  if(isCHadronFromB(parent))return true;
554  }
555  }
556 
557  return false;
558  }
559 
560  // Define the function fillHadronMap that fills the map of hadrons with their flavour.
561  void HadronOriginClassifier::fillHadronMap(std::set<const xAOD::TruthParticle*>& usedHadron, std::map<const xAOD::TruthParticle*,int>& mainHadronMap, const xAOD::TruthParticle* mainhad, const xAOD::TruthParticle* ihad, bool decayed) const {
562  // Fist, check that the consdired hadron has a non-null pointer
563  if (!ihad) return;
564  usedHadron.insert(ihad);
565  // Create two variables to indicate the flavour of the parents and childrens particles that will be considered.
566  // Create a boolean to indicate if the particles considered are from the final state.
567  int parent_flav,child_flav;
568  bool isFinal = true;
569  // Check if the considered hadron has children.
570  if(!ihad->nChildren()) return;
571  // Use a for to go through the children.
572  for(unsigned int j=0; j<ihad->nChildren(); ++j){
573  // Extract the j-th children.
574  const xAOD::TruthParticle* child = ihad->child(j);
575  if(!child) continue;
576  if(decayed){
577  fillHadronMap(usedHadron, mainHadronMap,mainhad,child,true);
578  isFinal=false;
579  }
580  else{
581  child_flav = std::abs(MC::leadingQuark(child));
582  if(child_flav!=4 && child_flav!=5) continue;
583  parent_flav = std::abs(MC::leadingQuark(mainhad));
584  if(child_flav!=parent_flav) continue;
585  fillHadronMap(usedHadron, mainHadronMap,mainhad,child);
586  isFinal=false;
587  }
588  }
589 
590  if(isFinal && !decayed){
591  mainHadronMap[mainhad]=std::abs(MC::leadingQuark(mainhad));
592  for(unsigned int j=0; j<ihad->nChildren(); ++j){
593  const xAOD::TruthParticle* child = ihad->child(j);
594  if(!child) continue;
595  fillHadronMap(usedHadron, mainHadronMap,mainhad,child,true);
596  }
597  }
598  }
599 
600  /*
601  ---------------------------------------------------------------------------------------------------------------------------------------
602  ----------------------------------------------------------- Particle Origin -----------------------------------------------------------
603  ---------------------------------------------------------------------------------------------------------------------------------------
604  */
605 
606  // Define the function isFromTop that indicates if a particle comes from top.
607 
609  // Find the first parent of the considered particle that is different from the particle.
610  const xAOD::TruthParticle* initpart = findInitial(part);
611  // Check if this parent comes from the top with function isDirectlyFromTop.
612  return isDirectlyFromTop(initpart);
613  }
614 
615  // Define the function isDirectlyFromTop that indicates if a particle comes from the direct decay of top.
617  // First, make sure the consdired particle has a non-null pointer and it has parents.
618  // Otherwise, return false.
619  if(!part || !part->nParents()) return false;
620  // Go through the parents of the particle.
621  for(unsigned int i=0; i<part->nParents(); ++i){
622  // Extract the i-th parent.
623  const xAOD::TruthParticle* parent = part->parent(i);
624  if(!parent) continue;
625  // If the i-th parent is a top, then return true
626  if( MC::isTop(parent) ) return true;
627  }
628  // If a top is no the parent, then return false.
629  return false;
630  }
631 
632  // Define the function isFromWTop that indicates if a particle comes from the decay chain t->Wb.
633 
635  // Find the first parent of the considered particle that is different from the particle.
636  const xAOD::TruthParticle* initpart = findInitial(part);
637  return isDirectlyFromWTop(initpart);
638  }
639 
640  // Define the function isDirectlyFromWTop that indicates if a particle comes from the direct decay of a W from a top.
642  // First, make sure the consdired particle has a non-null pointer and it has parents.
643  // Otherwise, return false.
644  if(!part || !part->nParents()) return false;
645  // Use a for to go though the parents.
646  for(unsigned int i=0; i<part->nParents(); ++i){
647  // Get the i-th parent.
648  const xAOD::TruthParticle* parent = part->parent(i);
649  if(!parent) continue;
650  if( MC::isW(parent)){
651  if( isFromTop(parent) ) return true;
652  }
653  }
654  // In this case, none of the parents of the particle is a W from top.
655  // Hence, return false.
656  return false;
657  }
658 
660  if(!part->nParents()) return false;
661  for(unsigned int i=0; i<part->nParents(); ++i){
662  const xAOD::TruthParticle* parent = part->parent(i);
663  if(!parent) continue;
664  if( MC::isPhoton(parent) || parent->absPdgId()<MC::BQUARK ) return true;
665  }
666  return false;
667  }
668 
670  const xAOD::TruthParticle* initpart = findInitial(part);
671  return isDirectlyFromGluonQuark(initpart);
672  }
673 
675  if(!part->nParents()) return false;
676  for(unsigned int i=0; i<part->nParents(); ++i){
677  const xAOD::TruthParticle* parent = part->parent(i);
678  if(!parent) continue;
679  if(!MC::isW(parent)) continue;
680  if(MC::isCharm(part)){
681  //trick to get at least 50% of PowhegPythia c from FSR
682  if(part->pdgId()==-(parent->pdgId())/6){
683  if( isFromGluonQuark(parent) ) return true;
684  }
685  }
686  else{
687  if( isFromGluonQuark(parent) ) return true;
688  }
689  }
690  return false;
691  }
692 
694  if(!part->nParents()) return false;
695  for(unsigned int i=0; i<part->nParents(); ++i){
696  const xAOD::TruthParticle* parent = part->parent(i);
697  if(!parent) continue;
699  if( isFromQuarkTop( parent ) ) return true;
700  }
701  }
702  return false;
703  }
704 
706  if(!part->nParents()) return false;
707  for(unsigned int i=0; i<part->nParents(); ++i){
708  const xAOD::TruthParticle* parent = part->parent(i);
709  if(!parent) continue;
710  if( parent->absPdgId()<MC::TQUARK ) {
711  if(isFromTop(parent)){
712  return true;
713  }
714  else if(isFromWTop(parent)){
715  return true;
716  }
717  }
718  }
719 
720  return false;
721  }
722 
724  const xAOD::TruthParticle* initpart = findInitial(part);
725  return isDirectlyFromQuarkTop(initpart);
726  }
727 
728  // Define the function isDirectlyFSRPythia8 that indicates if a particle comes from Final State Radiation in samples generated with Pythia8.
729 
731  // First, check if the particle has parents and return false if it does not.
732  if(!part->nParents()) return false;
733  // Use a for to go through the parents.
734  for(unsigned int i=0; i<part->nParents(); ++i){
735 
736  // Extract the i-th parent.
737 
738  const xAOD::TruthParticle* parent = part->parent(i);
739  if(!parent) continue;
741  if( isFromQuarkTopPythia8( parent ) ) return true;
742  }
743  }
744  // In this case, no parent from the particle is a gluon or a photon coming from a top
745  // Hence, the particle is not from FSR and false is not returned.
746  return false;
747  }
748 
749  // Define the function isDirectlyFromQuarkTopPythia8 that indicates if a particle comes from direct decay of the top in samples generated with Pythia8.
751  // First, make sure the consdired particle has a non-null pointer and it has parents.
752  // Otherwise, return false.
753  if(!part->nParents()) return false;
754  // Use a for to go through the parents.
755  for(unsigned int i=0; i<part->nParents(); ++i){
756  // Extract the i-th parent.
757  const xAOD::TruthParticle* parent = part->parent(i);
758  if(!parent) continue;
759  // Check if the parent is a quark different from the top.
760  if( parent->absPdgId()<MC::TQUARK ) {
761  // In this case, the parent is a quark different from top.
762  // Check if it comes from the decay chain of the t->Wb.
763  // If it is the case, return true.
764  if(isFromWTop(parent)){
765  return true;
766  }
767  }
768  }
769  // In this case, any of the parents of the particle comes from t->Wb chaing.
770  // Hence, the particle does not come from the top directly and false is returned.
771  return false;
772  }
773 
774  // Define the function isFromQuarkTopPythia8 that indicates if a particle comes from top in samples generated with Pythia8.
776  // Find the first parent of the considered particle that is different from the particle.
777  const xAOD::TruthParticle* initpart = findInitial(part);
778  // Check if this parent comes from the top with function isDirectlyFromQuarkTopPythia8.
779  return isDirectlyFromQuarkTopPythia8(initpart);
780  }
781 
783  if(!part->nParents()) return false;
784  for(unsigned int i=0; i<part->nParents(); ++i){
785  const xAOD::TruthParticle* parent = part->parent(i);
786  if(!parent) continue;
787  if( parent->absPdgId() == MC::PROTON && MC::isPhysical(part) ) return true;
788  }
789  return false;
790  }
791 
793  const xAOD::TruthParticle* initpart = findInitial(part);
794  return MC::Pythia8::isConditionC(initpart);
795  }
796 
798  if(!part->hasProdVtx()) return false;
799  const xAOD::TruthVertex* vertex = part->prodVtx();
800  return HepMC::status(vertex) == 2;
801  }
802 
803  /*
804  --------------------------------------------------------------------------------------------------------------------------------------
805  ---------------------------------------------------------- Particle Parents ----------------------------------------------------------
806  --------------------------------------------------------------------------------------------------------------------------------------
807  */
808 
809  // Define the function findInitial which finds the first parent of a particle that is not the particle itself.
810  const xAOD::TruthParticle* HadronOriginClassifier::findInitial(const xAOD::TruthParticle* part, std::shared_ptr<std::set<const xAOD::TruthParticle*>> checked ) const{
811  // If the particle has no parent, return the particle.
812  if(!part->nParents()) return part;
813  if (!checked) checked = std::make_shared<std::set<const xAOD::TruthParticle*>>();
814  // Use a for to go through the parents.
815  for(unsigned int i=0; i<part->nParents(); ++i){
816  // Extract the i-th parent.
817  const xAOD::TruthParticle* parent = part->parent(i);
818  if(!parent) continue;
819  if(checked->count(parent)) continue;
820  checked->insert(parent);
821  // If the parent has the same pdgId as the particle, then it means that the parent is the same as the considered particle.
822  // This happens if the particle irradiates for example.
823  // In this case, try to look for the first parent of i-th parent that is being considered.
824  // Repeat the process until you find a particle different from the considred one or that has no parent.
825 
826  if( part->pdgId() == parent->pdgId() ){
827  return findInitial(parent, std::move(checked));
828  }
829  }
830  // In this case, no parent different from the considered particle has been found.
831  // Hence, return the particle.
832  return part;
833  }
834 
835 }//namespace
LArG4FSStartPointFilter.part
part
Definition: LArG4FSStartPointFilter.py:21
xAOD::iterator
JetConstituentVector::iterator iterator
Definition: JetConstituentVector.cxx:68
DerivationFramework::HadronOriginClassifier::isDirectlyMPIPythia8
bool isDirectlyMPIPythia8(const xAOD::TruthParticle *part) const
Definition: HadronOriginClassifier.cxx:792
DerivationFramework::HadronOriginClassifier::m_DSID
Gaudi::Property< int > m_DSID
Definition: HadronOriginClassifier.h:90
MCTruthPartClassifier::hadron
@ hadron
Definition: TruthClassifiers.h:148
DerivationFramework::HadronOriginClassifier::isDirectlyMPISherpa
static bool isDirectlyMPISherpa(const xAOD::TruthParticle *part)
Definition: HadronOriginClassifier.cxx:797
ATH_MSG_INFO
#define ATH_MSG_INFO(x)
Definition: AthMsgStreamMacros.h:31
find
std::string find(const std::string &s)
return a remapped string
Definition: hcg.cxx:135
DerivationFramework::HadronOriginClassifier::GEN_id::HerwigPP
@ HerwigPP
SG::ReadHandle
Definition: StoreGate/StoreGate/ReadHandle.h:67
DerivationFramework::HadronOriginClassifier::b_FSR
@ b_FSR
Definition: HadronOriginClassifier.h:39
HadronOriginClassifier.h
DerivationFramework::HadronOriginClassifier::isFromQuarkTop
bool isFromQuarkTop(const xAOD::TruthParticle *part) const
Definition: HadronOriginClassifier.cxx:723
DerivationFramework::HadronOriginClassifier::GEN_id::Sherpa
@ Sherpa
MC::Pythia8::isConditionC
bool isConditionC(const T &p)
Definition: HepMCHelpers.h:27
DerivationFramework::HadronOriginClassifier::c_FSR
@ c_FSR
Definition: HadronOriginClassifier.h:39
DerivationFramework::HadronOriginClassifier::IsHerwigPP
bool IsHerwigPP() const
Definition: HadronOriginClassifier.h:82
skel.it
it
Definition: skel.GENtoEVGEN.py:407
DerivationFramework::HadronOriginClassifier::isFromQuarkTopPythia8
bool isFromQuarkTopPythia8(const xAOD::TruthParticle *part) const
Definition: HadronOriginClassifier.cxx:775
DerivationFramework::HadronOriginClassifier::GetOriginMap
std::map< const xAOD::TruthParticle *, HF_id > GetOriginMap() const
Definition: HadronOriginClassifier.cxx:360
DerivationFramework::HadronOriginClassifier::c_MPI
@ c_MPI
Definition: HadronOriginClassifier.h:38
read_hist_ntuple.t
t
Definition: read_hist_ntuple.py:5
DerivationFramework::HadronOriginClassifier::c_from_W
@ c_from_W
Definition: HadronOriginClassifier.h:40
DerivationFramework::HadronOriginClassifier::isFromTop
bool isFromTop(const xAOD::TruthParticle *part) const
Definition: HadronOriginClassifier.cxx:608
isBottomHadron
bool isBottomHadron(const T &p)
Definition: AtlasPID.h:906
master.gen
gen
Definition: master.py:32
DerivationFramework::HadronOriginClassifier::c_from_top
@ c_from_top
Definition: HadronOriginClassifier.h:41
isGluon
bool isGluon(const T &p)
Definition: AtlasPID.h:370
DerivationFramework::HadronOriginClassifier::m_mcName
SG::ReadHandleKey< xAOD::TruthEventContainer > m_mcName
Definition: HadronOriginClassifier.h:87
DerivationFramework::HadronOriginClassifier::isDirectlyFSR
bool isDirectlyFSR(const xAOD::TruthParticle *part) const
Definition: HadronOriginClassifier.cxx:693
DerivationFramework::HadronOriginClassifier::b_from_W
@ b_from_W
Definition: HadronOriginClassifier.h:40
DerivationFramework::HadronOriginClassifier::fillHadronMap
void fillHadronMap(std::set< const xAOD::TruthParticle * > &usedHadron, std::map< const xAOD::TruthParticle *, int > &mainHadronMap, const xAOD::TruthParticle *mainhad, const xAOD::TruthParticle *ihad, bool decayed=false) const
Definition: HadronOriginClassifier.cxx:561
DerivationFramework::HadronOriginClassifier::findInitial
const xAOD::TruthParticle * findInitial(const xAOD::TruthParticle *part, std::shared_ptr< std::set< const xAOD::TruthParticle * >> checked=nullptr) const
Definition: HadronOriginClassifier.cxx:810
DerivationFramework::HadronOriginClassifier::~HadronOriginClassifier
virtual ~HadronOriginClassifier()
Definition: HadronOriginClassifier.cxx:34
MC::isPhysical
bool isPhysical(const T &p)
Identify if the particle is physical, i.e. is stable or decayed.
Definition: HepMCHelpers.h:51
DerivationFramework::HadronOriginClassifier::IsSherpa
bool IsSherpa() const
Definition: HadronOriginClassifier.h:85
python.utils.AtlRunQueryDQUtils.p
p
Definition: AtlRunQueryDQUtils.py:209
isBottom
bool isBottom(const T &p)
Definition: AtlasPID.h:179
HepMC::is_simulation_particle
bool is_simulation_particle(const T &p)
Method to establish if a particle (or barcode) was created during the simulation (TODO update to be s...
Definition: MagicNumbers.h:354
DerivationFramework::HadronOriginClassifier::b_MPI
@ b_MPI
Definition: HadronOriginClassifier.h:38
lumiFormat.i
int i
Definition: lumiFormat.py:85
DerivationFramework::HadronOriginClassifier::isFromGluonQuark
bool isFromGluonQuark(const xAOD::TruthParticle *part) const
Definition: HadronOriginClassifier.cxx:669
leadingQuark
int leadingQuark(const T &p)
Definition: AtlasPID.h:883
DerivationFramework::HadronOriginClassifier::isDirectlyFromQuarkTop
bool isDirectlyFromQuarkTop(const xAOD::TruthParticle *part) const
Definition: HadronOriginClassifier.cxx:705
beamspotman.n
n
Definition: beamspotman.py:727
EL::StatusCode
::StatusCode StatusCode
StatusCode definition for legacy code.
Definition: PhysicsAnalysis/D3PDTools/EventLoop/EventLoop/StatusCode.h:22
DerivationFramework::HadronOriginClassifier::isFromWTop
bool isFromWTop(const xAOD::TruthParticle *part) const
Definition: HadronOriginClassifier.cxx:634
DerivationFramework::HadronOriginClassifier::GEN_id
GEN_id
Definition: HadronOriginClassifier.h:44
xAOD::TruthParticle_v1
Class describing a truth particle in the MC record.
Definition: TruthParticle_v1.h:37
DerivationFramework::HadronOriginClassifier::isDirectlyFromWTop
bool isDirectlyFromWTop(const xAOD::TruthParticle *part) const
Definition: HadronOriginClassifier.cxx:641
DerivationFramework::HadronOriginClassifier::isDirectlyFromGluonQuark
static bool isDirectlyFromGluonQuark(const xAOD::TruthParticle *part)
Definition: HadronOriginClassifier.cxx:659
DerivationFramework::HadronOriginClassifier::m_GenUsed
GEN_id m_GenUsed
Definition: HadronOriginClassifier.h:91
test_pyathena.parent
parent
Definition: test_pyathena.py:15
DerivationFramework::HadronOriginClassifier::GEN_id::Pythia8
@ Pythia8
DerivationFramework::HadronOriginClassifier::isDirectlyFSRPythia6
bool isDirectlyFSRPythia6(const xAOD::TruthParticle *part) const
Definition: HadronOriginClassifier.cxx:674
ATH_CHECK
#define ATH_CHECK
Definition: AthCheckMacros.h:40
DerivationFramework
THE reconstruction tool.
Definition: ParticleSortingAlg.h:24
DerivationFramework::HadronOriginClassifier::HadronOriginClassifier
HadronOriginClassifier(const std::string &t, const std::string &n, const IInterface *p)
Definition: HadronOriginClassifier.cxx:29
SG::ReadHandle::isValid
virtual bool isValid() override final
Can the handle be successfully dereferenced?
xAOD::TruthVertex_v1
Class describing a truth vertex in the MC record.
Definition: TruthVertex_v1.h:37
xAOD::TruthParticle_v1::nChildren
size_t nChildren() const
Number of children of this particle.
Definition: TruthParticle_v1.cxx:135
id
SG::auxid_t id
Definition: Control/AthContainers/Root/debug.cxx:239
DerivationFramework::HadronOriginClassifier::IsPythia6
bool IsPythia6() const
Definition: HadronOriginClassifier.h:84
xAOD::TruthParticle_v1::child
const TruthParticle_v1 * child(size_t i) const
Retrieve the i-th mother (TruthParticle) of this TruthParticle.
Definition: TruthParticle_v1.cxx:144
DerivationFramework::HadronOriginClassifier::isDirectlyFSRPythia8
bool isDirectlyFSRPythia8(const xAOD::TruthParticle *part) const
Definition: HadronOriginClassifier.cxx:730
DerivationFramework::HadronOriginClassifier::isDirectlyMPIPythia6
static bool isDirectlyMPIPythia6(const xAOD::TruthParticle *part)
Definition: HadronOriginClassifier.cxx:782
python.PyAthena.v
v
Definition: PyAthena.py:154
DerivationFramework::HadronOriginClassifier::IsPythia8
bool IsPythia8() const
Definition: HadronOriginClassifier.h:83
isW
bool isW(const T &p)
Definition: AtlasPID.h:379
Trk::vertex
@ vertex
Definition: MeasurementType.h:21
isTop
bool isTop(const T &p)
Definition: AtlasPID.h:182
ATH_MSG_WARNING
#define ATH_MSG_WARNING(x)
Definition: AthMsgStreamMacros.h:32
xAOD::EgammaHelpers::isPhoton
bool isPhoton(const xAOD::Egamma *eg)
is the object a photon
Definition: EgammaxAODHelpers.cxx:21
isCharm
bool isCharm(const T &p)
Definition: AtlasPID.h:176
DerivationFramework::HadronOriginClassifier::isDirectlyFromTop
static bool isDirectlyFromTop(const xAOD::TruthParticle *part)
Definition: HadronOriginClassifier.cxx:616
isCharmHadron
bool isCharmHadron(const T &p)
Definition: AtlasPID.h:905
python.SystemOfUnits.s
float s
Definition: SystemOfUnits.py:147
DerivationFramework::HadronOriginClassifier::isDirectlyFromQuarkTopPythia8
bool isDirectlyFromQuarkTopPythia8(const xAOD::TruthParticle *part) const
Definition: HadronOriginClassifier.cxx:750
DerivationFramework::HadronOriginClassifier::extrajet
@ extrajet
Definition: HadronOriginClassifier.h:37
ReadHandle.h
Handle class for reading from StoreGate.
AthAlgTool
Definition: AthAlgTool.h:26
HepMC::status
int status(const T &p)
Definition: MagicNumbers.h:143
DerivationFramework::HadronOriginClassifier::buildPartonsHadronsMaps
void buildPartonsHadronsMaps(std::map< const xAOD::TruthParticle *, int > &mainHadronMap, std::map< const xAOD::TruthParticle *, HF_id > &partonsOrigin) const
Definition: HadronOriginClassifier.cxx:428
HepMCHelpers.h
DerivationFramework::HadronOriginClassifier::b_from_top
@ b_from_top
Definition: HadronOriginClassifier.h:41
DerivationFramework::HadronOriginClassifier::isCHadronFromB
bool isCHadronFromB(const xAOD::TruthParticle *part, std::shared_ptr< std::set< const xAOD::TruthParticle * >> checked=nullptr) const
Definition: HadronOriginClassifier.cxx:539
DerivationFramework::HadronOriginClassifier::GEN_id::Pythia6
@ Pythia6
DerivationFramework::HadronOriginClassifier::initialize
virtual StatusCode initialize() override
Definition: HadronOriginClassifier.cxx:36