ATLAS Offline Software
Loading...
Searching...
No Matches
CalcPartonHistory.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
7
9
11
18
19#ifdef XAOD_STANDALONE
20#define TDS() evtStore()->tds()
21#else
22#define TDS() evtStore()
23#endif
24
25namespace {
26const std::vector<const xAOD::TruthParticle*>* findVector(
27 const std::map<std::string, std::vector<const xAOD::TruthParticle*>>& map,
28 const std::string& key) {
29 auto it = map.find(key);
30 return it == map.end() ? nullptr : &it->second;
31}
32
33constexpr const char* kBeforeFSR = "_beforeFSR";
34constexpr const char* kAfterFSR = "_afterFSR";
35} // namespace
36
37namespace CP {
38using ROOT::Math::PtEtaPhiMVector;
39
41 const std::string& name, const std::vector<std::string>& truthCollection)
42 : asg::AsgTool(name), m_truthCollections(truthCollection) {
43 declareProperty("prefix", m_prefix = "",
44 "Prefix to apply to all names to avoid overwriting");
45}
46
47bool CalcPartonHistory::ExistsInMap(const std::string& key) const {
48 // Checks whether a given key exists in the particle map.
49 return m_particleMap.find(key) != m_particleMap.end();
50}
51
52bool CalcPartonHistory::ExistsInKey(const std::string& key,
53 const xAOD::TruthParticle* p) const {
54 // Checks whether a given particle exists in the vector for the given key.
55 if (const auto* v = findVector(m_particleMap, key)) {
56 return std::find(v->begin(), v->end(), p) != v->end();
57 }
58 return false;
59}
60
61bool CalcPartonHistory::Retrievep4(const std::string& key,
62 PtEtaPhiMVector& p4) {
63 return Retrievep4(key, p4, 0);
64}
65
66bool CalcPartonHistory::Retrievep4(const std::string& key, PtEtaPhiMVector& p4,
67 const int& idx) {
68 const auto* v = findVector(m_particleMap, key);
69 if (!v || idx < 0)
70 return false;
71 const auto uidx = static_cast<std::size_t>(idx);
72 if (uidx >= v->size())
73 return false;
74 p4 = GetPtEtaPhiMfromTruth(v->at(uidx));
75 return true;
76}
77
78bool CalcPartonHistory::Retrievep4Gamma(PtEtaPhiMVector& p4, int& parentpdgId) {
79 // Finds the highest-pT photon from any 'GammaRad' key in the m_particleMap.
80 const xAOD::TruthParticle* bestPhoton = nullptr;
81 for (const auto& [key, particles] : m_particleMap) {
82 if (key.find("GammaRad") != std::string::npos) {
83 auto it = std::max_element(
84 particles.begin(), particles.end(),
85 [](const auto* a, const auto* b) { return a->pt() < b->pt(); });
86 if (it != particles.end() &&
87 (!bestPhoton || (*it)->pt() > bestPhoton->pt()))
88 bestPhoton = *it;
89 }
90 }
91 if (!bestPhoton) {
92 parentpdgId = 0;
93 return false;
94 }
95 p4 = GetPtEtaPhiMfromTruth(bestPhoton);
96 parentpdgId = bestPhoton->nParents() > 0 ? bestPhoton->parent(0)->pdgId() : 0;
97 return true;
98}
99
100bool CalcPartonHistory::RetrievepdgId(const std::string& key, int& pdgId) {
101 return RetrievepdgId(key, pdgId, 0);
102}
103
104bool CalcPartonHistory::RetrievepdgId(const std::string& key,
105 std::vector<int>& pdgIds) {
106 const auto* v = findVector(m_particleMap, key);
107 if (!v)
108 return false;
109 pdgIds.reserve(pdgIds.size() + v->size());
110 for (const auto* p : *v)
111 if (p)
112 pdgIds.push_back(p->pdgId());
113 return true;
114}
115
116bool CalcPartonHistory::RetrievepdgId(const std::string& key, int& pdgId,
117 const int& idx) {
118 const auto* v = findVector(m_particleMap, key);
119 if (!v || idx < 0)
120 return false;
121 const auto uidx = static_cast<std::size_t>(idx);
122 if (uidx >= v->size())
123 return false;
124 pdgId = v->at(uidx)->pdgId();
125 return true;
126}
127
129 const std::string& prefix, std::vector<const xAOD::TruthParticle*>& out) {
130 const auto* v = findVector(m_particleMap, prefix);
131 if (!v)
132 return false;
133 out.insert(out.end(), v->begin(), v->end());
134 return true;
135}
136
138 const std::string& prefix, std::vector<PtEtaPhiMVector>& particles,
139 std::vector<int>& pdgIds) {
140 const auto* v = findVector(m_particleMap, prefix);
141 if (!v)
142 return false;
143 particles.reserve(particles.size() + v->size());
144 pdgIds.reserve(pdgIds.size() + v->size());
145 for (const auto* p : *v) {
146 particles.push_back(GetPtEtaPhiMfromTruth(p));
147 pdgIds.push_back(p->pdgId());
148 }
149 return true;
150}
151
152bool CalcPartonHistory::RetrieveParticleInfo(const std::string& prefix,
153 PtEtaPhiMVector& particle,
154 int& pdgId) {
155 return RetrieveParticleInfo(prefix, particle, pdgId, 0);
156}
157
158bool CalcPartonHistory::RetrieveParticleInfo(const std::string& prefix,
159 PtEtaPhiMVector& particle,
160 int& pdgId, const int& idx) {
161 return Retrievep4(prefix, particle, idx) && RetrievepdgId(prefix, pdgId, idx);
162}
163
164bool CalcPartonHistory::RetrieveParticleInfo(const std::string& prefix,
165 const std::string& alt_prefix,
166 PtEtaPhiMVector& particle,
167 int& pdgId) {
168 if (Retrievep4(prefix, particle) && RetrievepdgId(prefix, pdgId))
169 return true;
170 return Retrievep4(alt_prefix, particle) && RetrievepdgId(alt_prefix, pdgId);
171}
172
174 const xAOD::TruthParticle* particle) {
175 static const std::unordered_map<int, std::string> pdgMap = {
176 {1, "_q"}, {2, "_q"}, {3, "_q"}, {-1, "_qbar"},
177 {-2, "_qbar"}, {-3, "_qbar"}, {6, "_t"}, {-6, "_tbar"},
178 {5, "_b"}, {-5, "_bbar"}, {4, "_c"}, {-4, "_cbar"},
179 {25, "_H"}, {24, "_W"}, {-24, "_W"}, {23, "_Z"},
180 {22, "_gamma"}, {21, "_g"}, {11, "_l"}, {13, "_l"},
181 {15, "_l"}, {-11, "_lbar"}, {-13, "_lbar"}, {-15, "_lbar"},
182 {12, "_nu"}, {14, "_nu"}, {16, "_nu"}, {-12, "_nubar"},
183 {-14, "_nubar"}, {-16, "_nubar"}, {2212, "_p"}, {1103, "_dd"},
184 {2101, "_ud"}, {2103, "_ud"}, {2203, "_uu"}, {3101, "_sd"},
185 {3103, "_sd"}, {3201, "_su"}, {3203, "_su"}, {3303, "_ss"},
186 {4101, "_cd"}, {4103, "_cd"}, {4201, "_cu"}, {4203, "_cu"},
187 {4301, "_cs"}, {4303, "_cs"}, {4403, "_cc"}, {5101, "_bd"},
188 {5103, "_bd"}, {5201, "_bu"}, {5203, "_bu"}, {5301, "_bs"},
189 {5303, "_bs"}, {5401, "_bc"}, {5403, "_bc"}, {5503, "_bb"}};
190 int pdgId = particle->pdgId();
191 auto it = pdgMap.find(pdgId);
192 return it != pdgMap.end() ? it->second : "_" + std::to_string(pdgId);
193}
194
196 const xAOD::TruthParticle* p, std::vector<const xAOD::TruthParticle*>& path,
197 std::vector<std::vector<const xAOD::TruthParticle*>>& allPaths) {
198 // Recursively builds decay paths from p down to stable particles.
199 // Each completed root-to-leaf path is appended to allPaths and later
200 // processed by FillParticleMap to assign m_particleMap keys.
201 //
202 // FSR handling: a particle that radiates before decaying appears as a chain
203 // of identical PDG-ID nodes (e.g. t → t → t → W b). findAfterFSR() follows
204 // that chain to the last node before the actual decay vertex, so the path
205 // records the pre- and post-FSR instances rather than every intermediate.
206 //
207 // W linking (DAOD_PHYS only): TruthTop links to Ws in TruthBoson, which
208 // carry no decay products. When we encounter such a W (isAfterFSR &&
209 // |pdg|==24) we swap it for the corresponding entry in
210 // TruthBosonsWithDecayParticles via the
211 // "CustomLinkedTruthBosonWithDecayParticles" decoration so that the subsequent
212 // child traversal finds the W decay products.
214 return;
215
216 // If this W node has no identical child (i.e. it is the after-FSR instance)
217 // but lacks decay products, swap to the linked decay-product-bearing copy.
218 if (PartonHistoryUtils::isAfterFSR(p) && std::abs(p->pdgId()) == 24) {
220 p, "CustomLinkedTruthBosonWithDecayParticles");
221 }
222
223 path.push_back(p);
224
225 // Leaf node: record the completed path.
226 if (p->nChildren() == 0) {
227 allPaths.push_back(path);
228 path.pop_back();
229 return;
230 }
231
233 // Same W-linking fix for the after-FSR node reached by findAfterFSR().
234 if (std::abs(afterFSR->pdgId()) == 24) {
236 afterFSR, "CustomLinkedTruthBosonWithDecayParticles");
237 }
238
239 if (afterFSR != p) {
240 // There was FSR: continue the path from the post-FSR node.
241 TraceParticle(afterFSR, path, allPaths);
242 } else {
243 // No FSR: branch into each child (e.g. W→lν gives two sub-paths).
244 for (std::size_t i = 0; i < afterFSR->nChildren(); ++i) {
245 if (const auto* c = afterFSR->child(i))
246 TraceParticle(c, path, allPaths);
247 }
248 }
249 path.pop_back();
250}
251
253 const std::string& key) {
254 if (!ExistsInKey(key, p))
255 m_particleMap[key].push_back(p);
256}
257
259 const std::string& newKey, std::string& key) {
261 key += newKey;
262
263 if (p->nParents() == 0) {
264 AddToParticleMap(p, key + kBeforeFSR);
266 AddToParticleMap(p, key + kAfterFSR);
267 return true;
268 }
269
271 AddToParticleMap(p, key + kAfterFSR);
272 } else {
273 AddToParticleMap(p, key + kBeforeFSR);
275 AddToParticleMap(p, key + kAfterFSR);
276 }
277 return true;
278}
279
281 std::string& key, int decayID) {
282 const bool fromH = PartonHistoryUtils::hasParentAbsPdgId(p, 25) &&
284 const bool fromW = PartonHistoryUtils::hasParentAbsPdgId(p, 24) &&
286 const bool fromZ = PartonHistoryUtils::hasParentAbsPdgId(p, 23) &&
288 if (!fromH && !fromW && !fromZ)
289 return false;
290
291 const std::string decayStr = "Decay" + std::to_string(decayID);
292 key += decayStr;
293 AddToParticleMap(p, key + kBeforeFSR);
295 AddToParticleMap(p, key + kAfterFSR);
296 return true;
297}
298
300 std::string& key) {
301 AddToParticleMap(particle, key);
302}
303
305 const std::string& newKey,
306 std::string& key) {
307 AddToParticleMap(particle, key + newKey);
308 key += newKey;
309}
310
312 std::vector<std::vector<const xAOD::TruthParticle*>>& allPaths) {
313 // Converts the raw decay paths produced by TraceParticles into the
314 // m_particleMap used by all Retrieve* and Fill* methods.
315 //
316 // Key construction: each path is walked particle by particle, accumulating
317 // a string key of the form
318 // "<prefix>_MC_<type1>_<type2>_..._<beforeFSR|afterFSR>". For example, a b
319 // quark from a top gives "MySch_MC_t_b_beforeFSR". Decay products of W/Z/H
320 // get an additional "Decay<N>" segment to distinguish the two daughters, e.g.
321 // "MySch_MC_t_WDecay1_beforeFSR".
322 //
323 // Handler priority (first match wins for each particle in the path):
324 // handleDecay — daughters of W/Z/H: appends "Decay<N>" and records
325 // beforeFSR/afterFSR handleFSR — the radiating particle itself:
326 // records both before and after FSR copies handleSameAsParent —
327 // intermediate FSR copies (same PDG as parent): stored under current key
328 // handleDefault — all other particles: appends the type suffix and
329 // advances the key
330 m_particleMap.clear();
331 static const SG::Accessor<unsigned int> acc_classification("Classification");
332 static const SG::Accessor<unsigned int> acc_classifierParticleOrigin(
333 "classifierParticleOrigin");
334 static const SG::Accessor<unsigned int> acc_classifierParticleType(
335 "classifierParticleType");
336
337 for (const auto& path : allPaths) {
338 // m_particleMap keys always include the prefix, built once here.
339 std::string key = m_prefix + "_" + "MC";
340
341 for (const auto* p : path) {
342 // beforeFSR: this node has an identical child (it will radiate).
343 // afterFSR: this node's parent has the same PDG ID (it was radiated
344 // from).
345 const bool beforeFSR = PartonHistoryUtils::hasIdenticalChild(p);
346 const bool afterFSR = PartonHistoryUtils::hasParentPdgId(p);
347
348 // Determine which child index this particle is under its parent; used to
349 // label W/Z/H decay daughters as Decay1, Decay2. Falls back to sign of
350 // pdgId (negative → 2) if parent navigation is unavailable.
351 int decayID = (p->pdgId() < 0) ? 2 : 1;
352 if (p->nParents() != 0 && p->parent(0)) {
353 const auto* par = p->parent(0);
354 for (std::size_t i = 0; i < par->nChildren(); ++i) {
355 if (par->child(i) == p) {
356 decayID = static_cast<int>(i) + 1;
357 break;
358 }
359 }
360 }
361
362 const std::string new_key = GetParticleType(p);
363
364 // Skip pure intermediate FSR nodes that are neither the first nor last
365 // in the FSR chain — they carry no additional physics information.
366 if (beforeFSR && afterFSR)
367 continue;
368 if (handleDecay(p, key, decayID))
369 continue;
370 if (handleFSR(p, new_key, key))
371 continue;
373 handleSameAsParent(p, key);
374 continue;
375 }
376 if (!new_key.empty())
377 handleDefault(p, new_key, key);
378 }
379 }
380}
381
383 const xAOD::TruthParticleContainer* truthParticles) {
384 std::vector<std::vector<const xAOD::TruthParticle*>> allPaths;
385 allPaths.reserve(truthParticles->size());
386 for (const xAOD::TruthParticle* p : *truthParticles) {
388 continue;
389 std::vector<const xAOD::TruthParticle*> path;
390 path.reserve(16);
391 TraceParticle(p, path, allPaths);
392 }
393 FillParticleMap(allPaths);
394}
395
397 m_dec.setPrefix(m_prefix);
399 return StatusCode::SUCCESS;
400}
401
403 const xAOD::TruthParticleContainer* truthParticles{nullptr};
404 ANA_CHECK(linkTruthContainers(truthParticles));
405
406 const xAOD::EventInfo* partonHistory = nullptr;
407 ANA_CHECK(evtStore()->retrieve(partonHistory, "EventInfo"));
408
409 ANA_CHECK(runHistorySaver(truthParticles, partonHistory));
410 return StatusCode::SUCCESS;
411}
412
414 const std::vector<std::string>& collections,
415 const std::string& out_contName) {
416 // In DAOD_PHYS there is no single TruthParticles collection, so we merge
417 // several dedicated collections (e.g. TruthTop,
418 // TruthBosonsWithDecayParticles) into one working container stored in the
419 // event store.
420 //
421 // The merge must avoid double-counting: the same particle can appear in
422 // multiple collections (e.g. the b from a top appears in both TruthTop and
423 // TruthBottom). We therefore keep only "root" particles — those that are not
424 // a descendant of any other candidate in the merged pool. TraceParticles then
425 // walks down from these roots, naturally visiting all descendants regardless
426 // of which original collection they came from.
430 std::vector<const xAOD::TruthParticle*> p_candidates;
431 std::vector<const xAOD::TruthParticle*> p_parents;
432
433 for (const std::string& collection : collections) {
434 const xAOD::TruthParticleContainer* cont = nullptr;
435 ANA_CHECK(evtStore()->retrieve(cont, collection));
436 p_candidates.insert(p_candidates.end(), cont->begin(), cont->end());
437 }
438 // Retain only particles that have no ancestor among the other candidates.
439 for (const xAOD::TruthParticle* potential_parent : p_candidates) {
440 if (std::none_of(p_candidates.begin(), p_candidates.end(),
441 [&](const xAOD::TruthParticle* other_candidate) {
442 return other_candidate != potential_parent &&
443 PartonHistoryUtils::isChildOf(other_candidate,
444 potential_parent);
445 }))
446 p_parents.push_back(potential_parent);
447 }
448 out_cont->insert(out_cont->end(), p_parents.begin(), p_parents.end());
449 StatusCode save = TDS()->record(out_cont, out_contName);
450 if (!save)
451 return StatusCode::FAILURE;
452 return StatusCode::SUCCESS;
453}
454
457 "TruthBoson", "TruthBosonsWithDecayParticles",
458 "CustomLinkedTruthBosonWithDecayParticles");
459}
460
462 const std::string& collectionToDecorate,
463 const std::string& collectionToLink, const std::string& nameOfDecoration) {
464 const SG::Decorator<const xAOD::TruthParticle*> dec(nameOfDecoration);
465 const xAOD::TruthParticleContainer* cont1 = nullptr;
466 const xAOD::TruthParticleContainer* cont2 = nullptr;
467 ANA_CHECK(evtStore()->retrieve(cont1, collectionToDecorate));
468 ANA_CHECK(evtStore()->retrieve(cont2, collectionToLink));
469 for (const auto* p : *cont1) {
470 const xAOD::TruthParticle* link = nullptr;
471 for (const auto* q : *cont2) {
472 if (p->pdgId() == q->pdgId() && p->uid() == q->uid()) {
473 link = q;
474 break;
475 }
476 }
477 dec(*p) = link;
478 }
479 return StatusCode::SUCCESS;
480}
481
484 const xAOD::TruthParticle* part, const std::string& decorationName) {
486 if (!acc.isAvailable(*part))
487 return part;
488 const xAOD::TruthParticle* link = acc(*part);
489 return link ? link : part;
490}
491
493 const xAOD::TruthParticleContainer*& tp) {
494 const std::string key = m_prefix + "_TruthParticles";
496 const auto& collections =
497 m_configured ? m_config.truthCollections : m_truthCollections;
499 ANA_CHECK(evtStore()->retrieve(tp, key));
501 return StatusCode::SUCCESS;
502 }
503 ANA_CHECK(evtStore()->retrieve(tp, key));
504 return StatusCode::SUCCESS;
505}
506
508 m_config = config;
509 m_configured = true;
510}
511
513 if (!m_configured)
514 return;
515
516 for (const auto& group : m_config.decoratorGroups) {
517 switch (group) {
520 break;
523 break;
526 break;
529 break;
532 break;
535 break;
538 break;
541 break;
544 break;
547 break;
550 break;
553 break;
556 break;
559 break;
560 }
561 }
562
563 for (const auto& zw : m_config.decoratorZWs) {
564 if (zw.type == DecoratorZW::Z)
565 InitializeZDecorators(zw.count, zw.extended);
566 else
567 InitializeWDecorators(zw.count);
568 }
569
570 for (const auto& fill : m_config.genericFills) {
571 if (fill.isVector) {
572 m_dec.initializeVectorPtEtaPhiMDecorator(fill.decorationKey);
573 m_dec.initializeVectorIntDecorator(fill.decorationKey + "_pdgId");
574 } else {
575 m_dec.initializePtEtaPhiMDecorator(fill.decorationKey);
576 m_dec.initializeIntDecorator(fill.decorationKey + "_pdgId");
577 }
578 }
579}
580
582 const xAOD::TruthParticleContainer* truthParticles,
583 const xAOD::EventInfo* partonHistory) {
584
585 // Register the EventInfo for this event; all m_dec.decorate*() calls
586 // will write to it automatically.
587 m_dec.setEventInfo(partonHistory);
588
589 TraceParticles(truthParticles);
590
591 for (const auto& op : m_config.specialFills) {
592 switch (op.type) {
595 break;
598 break;
601 break;
603 FillZPartonHistory(op.parent, op.count, op.mode);
604 break;
606 FillZtautauPartonHistory(op.parent, op.count, op.mode);
607 break;
609 FillWPartonHistory(op.parent, op.count, op.mode);
610 break;
612 FillHiggsPartonHistory(op.mode);
613 break;
615 FillGammaPartonHistory(op.parent);
616 break;
617 }
618 }
619
620 for (const auto& fill : m_config.genericFills) {
621 if (fill.isVector) {
622 FillGenericVectorPartonHistory(fill.retrievalKeys.at(0),
623 fill.decorationKey);
624 } else if (fill.retrievalKeys.size() == 1) {
625 FillGenericPartonHistory(fill.retrievalKeys.at(0), fill.decorationKey,
626 fill.idx);
627 } else {
628 FillGenericPartonHistory(fill.retrievalKeys, fill.decorationKey,
629 fill.idx);
630 }
631 }
632
633 return StatusCode::SUCCESS;
634}
635
636} // namespace CP
#define TDS()
DataVector adapter that acts like it holds const pointers.
#define ANA_CHECK(EXP)
check whether the given expression was successful
static Double_t a
ROOT::Math::PtEtaPhiMVector GetPtEtaPhiMfromTruth(const xAOD::TruthParticle *TruthParticle)
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
ServiceHandle< StoreGateSvc > & evtStore()
StatusCode buildContainerFromMultipleCollections(const std::vector< std::string > &collections, const std::string &out_contName)
used to build container from multiple collections in DAOD_PHYS we don't have the TruthParticles colle...
CalcPartonHistory(const std::string &name, const std::vector< std::string > &truthCollections={"TruthTop"})
void FillGenericVectorPartonHistory(const std::string &retrievalstring, const std::string &decorationstring)
bool ExistsInKey(const std::string &key, const xAOD::TruthParticle *particle) const
void configure(const PartonSchemeConfig &config)
const xAOD::TruthParticle * getTruthParticleLinkedFromDecoration(const xAOD::TruthParticle *part, const std::string &decorationName)
helper method to handle retriveing the truth particle linked in the decoration of another particle
virtual StatusCode runHistorySaver(const xAOD::TruthParticleContainer *truthParticles, const xAOD::EventInfo *ttbarPartonHistory)
void FillParticleMap(std::vector< std::vector< const xAOD::TruthParticle * > > &allPaths)
bool Retrievep4Gamma(PtEtaPhiMVector &p4, int &parentpdgId)
bool RetrieveParticleInfo(const std::string &prefix, std::vector< const xAOD::TruthParticle * > &particles)
const std::vector< std::string > m_truthCollections
void FillHiggsPartonHistory(const std::string &mode)
std::map< std::string, std::vector< const xAOD::TruthParticle * > > m_particleMap
virtual StatusCode linkTruthContainers(const xAOD::TruthParticleContainer *&tp)
void FillWPartonHistory(const std::string &parent, int nWs=1, const std::string &mode="resonant")
std::string GetParticleType(const xAOD::TruthParticle *particle)
bool Retrievep4(const std::string &key, PtEtaPhiMVector &p4)
std::string m_prefix
prefix applied to all decorator and m_particleMap names
void FillGammaPartonHistory(const std::string &parent)
bool ExistsInMap(const std::string &key) const
void InitializeZDecorators(int nZs=1, bool extend=false)
void FillZtautauPartonHistory(const std::string &parent, int nZs=1, const std::string &mode="resonant")
bool RetrievepdgId(const std::string &key, std::vector< int > &pdgIds)
bool handleDecay(const xAOD::TruthParticle *particle, std::string &key, int decayID)
bool m_configured
true after configure() has been called
void FillZPartonHistory(const std::string &parent, int nZs=1, const std::string &mode="resonant")
void TraceParticle(const xAOD::TruthParticle *particle, std::vector< const xAOD::TruthParticle * > &currentPath, std::vector< std::vector< const xAOD::TruthParticle * > > &allPaths)
void handleDefault(const xAOD::TruthParticle *particle, const std::string &newKey, std::string &key)
virtual void initializeDecorators()
virtual StatusCode execute()
void TraceParticles(const xAOD::TruthParticleContainer *truthParticles)
StatusCode linkBosonCollections()
currently in DAOD_PHYS TruthTop have links to Ws from the TruthBoson collection, which have no link t...
void FillGenericPartonHistory(const std::string &retrievalstring, const std::string &decorationstring, const int idx)
StatusCode decorateCollectionWithLinksToAnotherCollection(const std::string &collectionToDecorate, const std::string &collectionToLink, const std::string &nameOfDecoration)
helper method currently used in DAOD_PHYS to link particles from a given collection to the same parti...
bool handleFSR(const xAOD::TruthParticle *particle, const std::string &newKey, std::string &key)
void handleSameAsParent(const xAOD::TruthParticle *particle, std::string &key)
PartonSchemeConfig m_config
scheme configuration set via configure()
void AddToParticleMap(const xAOD::TruthParticle *particle, const std::string &key)
virtual StatusCode initialize() override
Dummy implementation of the initialisation function.
DataVector adapter that acts like it holds const pointers.
iterator begin() noexcept
Return an iterator pointing at the beginning of the collection.
iterator end() noexcept
Return an iterator pointing past the end of the collection.
iterator insert(iterator position, value_type pElem)
Add a new element to the collection.
const_iterator end() const noexcept
Return a const_iterator pointing past the end of the collection.
const_iterator begin() const noexcept
Return a const_iterator pointing at the beginning of the collection.
size_type size() const noexcept
Returns the number of elements in the collection.
Helper class to provide type-safe access to aux data.
Helper class to provide constant type-safe access to aux data.
Helper class to provide type-safe access to aux data.
Definition Decorator.h:59
AsgTool(const std::string &name)
Constructor specifying the tool instance's name.
Definition AsgTool.cxx:58
STL class.
const TruthParticle_v1 * child(size_t i) const
Retrieve the i-th mother (TruthParticle) of this TruthParticle.
int pdgId() const
PDG ID code.
const TruthParticle_v1 * parent(size_t i) const
Retrieve the i-th mother (TruthParticle) of this TruthParticle.
virtual double pt() const override final
The transverse momentum ( ) of the particle.
size_t nParents() const
Number of parents of this particle.
size_t nChildren() const
Number of children of this particle.
bool contains(const std::string &s, const std::string &regx)
does a string contain the substring
Definition hcg.cxx:114
bool hasIdenticalChild(const xAOD::TruthParticle *particle)
bool hasParticleIdenticalParent(const xAOD::TruthParticle *particle)
Return true when particle is a top before FSR.
bool hasParentPdgId(const xAOD::TruthParticle *particle, int PdgId)
const xAOD::TruthParticle * findAfterFSR(const xAOD::TruthParticle *particle)
Return particle after FSR (before the decay vertex)
bool isAfterFSR(const xAOD::TruthParticle *particle)
Determine whether particle is afterFSR.
bool isBrokenTop(const xAOD::TruthParticle *particle)
Looking for tops without children -> must be broken.
bool hasParentAbsPdgId(const xAOD::TruthParticle *particle, int absPdgId)
Select isolated Photons, Electrons and Muons.
@ Z
FillZPartonHistory(history, parent, dec, count, mode)
@ Ttbar
FillTtbarPartonHistory.
@ W
FillWPartonHistory(history, parent, dec, count, mode)
@ Ztautau
FillZtautauPartonHistory(history, parent, dec, count, mode)
@ Higgs
FillHiggsPartonHistory(history, mode, dec)
@ Top
FillTopPartonHistory.
@ AntiTop
FillAntiTopPartonHistory.
@ Gamma
FillGammaPartonHistory(history, parent, dec)
@ FourTop
Initialize4TopDecorators()
@ VectorAntiBottom
InitializeVectorAntiBottomDecorators()
@ VectorAntiCharm
InitializeVectorAntiCharmDecorators()
@ Bottom
InitializeBottomDecorators()
@ Charm
InitializeCharmDecorators()
@ AntiCharm
InitializeAntiCharmDecorators()
@ VectorBottom
InitializeVectorBottomDecorators()
@ Ttbar
InitializeTtbarDecorators()
@ VectorCharm
InitializeVectorCharmDecorators()
@ Higgs
InitializeHiggsDecorators()
@ Top
InitializeTopDecorators()
@ AntiBottom
InitializeAntiBottomDecorators()
@ AntiTop
InitializeAntiTopDecorators()
@ Photon
InitializePhotonDecorators()
@ VIEW_ELEMENTS
this data object is a view, it does not own its elmts
EventInfo_v1 EventInfo
Definition of the latest event info version.
TruthParticle_v1 TruthParticle
Typedef to implementation.
static const SG::AuxElement::Accessor< ElementLink< IParticleContainer > > acc("originalObjectLink")
Object used for setting/getting the dynamic decoration in question.
TruthParticleContainer_v1 TruthParticleContainer
Declare the latest version of the truth particle container.
Top-level configuration for a named parton history scheme.
void fill(H5::Group &out_file, size_t iterations)