ATLAS Offline Software
Loading...
Searching...
No Matches
HepMCHelpers.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4#ifndef TRUTHUTILS_HEPMCHELPERS_H
5#define TRUTHUTILS_HEPMCHELPERS_H
6#include <vector>
7#include <cmath>
8#include <algorithm>
9#include <array>
10#include <cstdlib>
11#include <set>
12#include <memory>
14
17/// ATLAS-specific HepMC functions
18
19#if !defined(HEPMC3) && !defined(XAOD_ANALYSIS)
21#include <ranges>
22namespace MC {
23inline
24auto particles_in (const HepMC::GenVertex* p) {
25 return std::ranges::subrange (p->particles_in_const_begin(),
26 p->particles_in_const_end());
29#endif
31namespace MC
33 template <class VTX>
34 auto particles_in (const VTX* p) { return p->particles_in(); }
35 template <class VTX>
36 auto particles_in (const std::shared_ptr<VTX>& p) { return p->particles_in(); }
37
38 namespace Pythia8
39 {
41 template <class T> inline bool isConditionA(const T& p) { return p->status() == 62 || p->status() == 52 || p->status() == 21 || p->status() == 22;}
42
43 template <class T> inline bool isConditionB(const T& p) { return p->status() == 23;}
44
45 template <class T> inline bool isConditionC(const T& p) { return p->status() > 30 && p->status() < 40;}
46 }
47
48#include "AtlasPID.h"
49
51 template <class T> inline bool isInteracting(const T& p) { return isStrongInteracting<T>(p) || isEMInteracting<T>(p) || isGeantino<T>(p); }
52
54 template <class T> inline bool isChargedNonShowering(const T& p) { return (isMuon<T>(p) || isSUSY<T>(p)); }
55
57 template <class T> inline bool isBeam(const T& p) { return HepMC::status(p)%HepMC::SIM_STATUS_THRESHOLD == 4;}
58
60 template <class T> inline bool isDecayed(const T& p) { return HepMC::status(p)%HepMC::SIM_STATUS_THRESHOLD == 2;}
61
63 template <class T> inline bool isStable(const T& p) { return HepMC::status(p)%HepMC::SIM_STATUS_THRESHOLD == 1;}
64
66 template <class T> inline bool isFinalState(const T& p) { return HepMC::status(p)%HepMC::SIM_STATUS_THRESHOLD == 1 && !p->end_vertex();}
69 template <class T> inline bool isPhysical(const T& p) { return isStable<T>(p) || isDecayed<T>(p); }
72 template <class T> inline bool isGenStable(const T& p) { return isStable<T>(p) && !HepMC::is_simulation_particle<T>(p);}
75 template <class T> inline bool isSimStable(const T& p) { return isStable<T>(p) && !p->end_vertex() && HepMC::is_simulation_particle<T>(p);}
78 template <class T> inline bool isSimInteracting(const T& p) { return isGenStable<T>(p) && isInteracting<T>(p);}
83 template <class T> inline bool isStableOrSimDecayed(const T& p) {
84 const auto vertex = p->end_vertex();
85 return ( isStable<T>(p) || (isDecayed<T>(p) && (!vertex || HepMC::is_simulation_vertex(vertex))));
86 }
89 template <class T> inline bool isZeroEnergyPhoton(const T& p) { return isPhoton<T>(p) && p->e() == 0;}
92 template <class T> inline bool isSpecialNonInteracting(const T& p) {
93 const int apid = std::abs(p->pdg_id());
94 if (apid == NU_E || apid == NU_MU || apid == NU_TAU) return true; //< neutrinos
95 if (apid == 1000022 || apid == 1000024 || apid == 5100022) return true; // SUSY & KK photon and Z partners
96 if (apid == GRAVITON || apid == 1000039 || apid == 5000039) return true; //< gravitons: standard, SUSY and KK
97 if (apid == 9000001 || apid == 9000002 || apid == 9000003 || apid == 9000004 || apid == 9000005 || apid == 9000006) return true; //< exotic particles from monotop model
98 return false;
99 }
102
103 template <class T> T findMother(T thePart) {
104 auto partOriVert = thePart->production_vertex();
105 if (!partOriVert) return nullptr;
106
107 long partPDG = thePart->pdg_id();
108 long MotherPDG(0);
109
110 auto MothOriVert = partOriVert;
111 MothOriVert = nullptr;
112 T theMoth(nullptr);
114 size_t itr = 0;
115 do {
116 if (itr != 0) partOriVert = MothOriVert;
117 for ( const auto& p : particles_in(partOriVert) ) {
118 theMoth = p;
119 if (!theMoth) continue;
120 MotherPDG = theMoth->pdg_id();
121 MothOriVert = theMoth->production_vertex();
122 if (MotherPDG == partPDG) break;
124 itr++;
125 if (itr > 100) {
126 break;
128 } while (MothOriVert != nullptr && MotherPDG == partPDG && !HepMC::is_simulation_particle(thePart) &&
129 MothOriVert != partOriVert);
130 return theMoth;
131 }
132
134
135 template <class C, class T> T findMatching(C TruthContainer, T p) {
136 T ptrPart = nullptr;
137 if (!p) return ptrPart;
138 if constexpr (std::is_pointer_v<C> || HepMC::is_smart_ptr_v<C>){ //C is ptr
139 for (T truthParticle : *TruthContainer) {
140 if (HepMC::is_sim_descendant(p,truthParticle)) {
141 ptrPart = truthParticle;
142 break;
146 else {
147 for (T truthParticle : TruthContainer) {
148 if (HepMC::is_sim_descendant(p,truthParticle)) {
149 ptrPart = truthParticle;
150 break;
153 }
154 return ptrPart;
155 }
156
157
158 template <class T> void findParticleAncestors(T thePart, std::set<T>& allancestors) {
159 auto prodVtx = thePart->production_vertex();
160 if (!prodVtx) return;
161 for (const auto& theMother: prodVtx->particles_in()) {
162 if (!theMother) continue;
163 allancestors.insert(theMother);
164 findParticleAncestors(theMother, allancestors);
165 }
166 }
167
170 template <class T> void findParticleStableDescendants(T thePart, std::set<T>& allstabledescendants) {
171 auto endVtx = thePart->end_vertex();
172 if (!endVtx) return;
173 for (const auto& theDaughter: endVtx->particles_out()) {
174 if (!theDaughter) continue;
175 if (isStable(theDaughter) && !HepMC::is_simulation_particle(theDaughter)) {
176 allstabledescendants.insert(theDaughter);
178 findParticleStableDescendants(theDaughter, allstabledescendants);
179 }
185
186 template <class T> bool isHardScatteringVertex(T pVert) {
187 if (pVert == nullptr) return false;
188 T pV = pVert;
189 int numOfPartIn(0);
190 int pdg(0);
192 do {
193 pVert = pV;
194 auto incoming = pVert->particles_in();
195 numOfPartIn = incoming.size();
196 pdg = numOfPartIn && incoming.front() != nullptr ? incoming.front()->pdg_id() : 0;
197 pV = numOfPartIn && incoming.front() != nullptr ? incoming.front()->production_vertex() : nullptr;
198
199 } while (numOfPartIn == 1 && (std::abs(pdg) < 81 || std::abs(pdg) > 100) && pV != nullptr);
201 if (numOfPartIn == 2) {
202 auto incoming = pVert->particles_in();
203 if (incoming.at(0) && incoming.at(1) && (std::abs(incoming.at(0)->pdg_id()) < 7 || incoming.at(0)->pdg_id() == 21) && (std::abs(incoming.at(1)->pdg_id()) < 7 || incoming.at(1)->pdg_id() == 21)) return true;
205 return false;
211
212 template <class T, class U>
213 bool isFromHadron(T p, U hadron, bool &fromTau, bool &fromBSM) {
214 if (isHadron(p)&&!isBeam(p)) return true; // trivial case
215 auto vtx = p->production_vertex();
216 if (!vtx) return false;
217 bool fromHad = false;
218 for ( const auto& parent : particles_in(vtx) ) {
219 if (!parent) continue;
220 // should this really go into parton-level territory?
221 // probably depends where BSM particles are being decayed
222 fromBSM |= isBSM(parent);
223 if (!isPhysical(parent)) return false;
224 fromTau |= isTau(parent);
225 if (isHadron(parent)&&!isBeam(parent)) {
226 if (!hadron) hadron = parent; // assumes linear hadron parentage
227 return true;
229 fromHad |= isFromHadron(parent, hadron, fromTau, fromBSM);
230 }
231 return fromHad;
232 }
233
235
236
237 template <class T> auto findSimulatedEndVertex(T thePart) -> decltype(thePart->end_vertex()) {
238 decltype(thePart->end_vertex()) EndVert = thePart->end_vertex();
239 decltype(thePart->end_vertex()) pVert(nullptr);
240 if (EndVert != nullptr) {
241 do {
242 bool samePart = false;
243 pVert = nullptr;
244 auto outgoing = EndVert->particles_out();
245 auto incoming = EndVert->particles_in();
246 for (const auto& itrDaug: outgoing) {
247 if (!itrDaug) continue;
248 if ((( HepMC::is_same_generator_particle(itrDaug,thePart)) ||
249 // brem on generator level for tau
250 (outgoing.size() == 1 && incoming.size() == 1 &&
252 itrDaug->pdg_id() == thePart->pdg_id()) {
253 samePart = true;
254 pVert = itrDaug->end_vertex();
255 }
256 }
257 if (samePart) EndVert = pVert;
258 } while (pVert != nullptr && pVert != EndVert); // pVert!=EndVert to prevent Sherpa loop
259 }
260 return EndVert;
261 }
262
264
265 template <class V> auto findFinalStateParticles(V theVert) -> decltype(theVert->particles_out()) {
266 if (!theVert) return {};
267 decltype(theVert->particles_out()) finalStatePart;
268 auto outgoing = theVert->particles_out();
269 for (const auto& thePart: outgoing) {
270 if (!thePart) continue;
271 finalStatePart.push_back(thePart);
272 if (isStable(thePart)) continue;
273 V pVert = findSimulatedEndVertex(thePart);
274 if (pVert == theVert) break; // to prevent Sherpa loop
275 if (pVert != nullptr) {
276 auto vecPart = findFinalStateParticles<V>(pVert);
277 finalStatePart.insert(finalStatePart.end(),vecPart.begin(),vecPart.end());
278 }
280 return finalStatePart;
281 }
282
284#endif
struct color C
bool is_same_generator_particle(const T1 &p1, const T2 &p2)
Method to establish if two particles in the GenEvent actually represent the same generated particle.
bool is_simulation_vertex(const T &v)
Method to establish if the vertex was created during simulation (TODO migrate to be based on status).
constexpr int SIM_STATUS_THRESHOLD
Constant definiting the status threshold for simulated particles, eg. can be used to separate generat...
int status(const T &p)
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...
constexpr bool is_smart_ptr_v
bool is_sim_descendant(const T1 &p1, const T2 &p2)
Method to check if the first particle is a descendant of the second in the simulation,...
bool isConditionB(const T &p)
bool isConditionA(const T &p)
To be understood.
bool isConditionC(const T &p)
T findMatching(C TruthContainer, T p)
Function to find a particle in container.
static const int GRAVITON
bool isZeroEnergyPhoton(const T &p)
Identify a photon with zero energy. Probably a workaround for a generator bug.
bool isHardScatteringVertex(T pVert)
Function to classify the vertex as hard scattering vertex.
bool isStableOrSimDecayed(const T &p)
Identify if particle is satble or decayed in simulation.
bool isPhoton(const T &p)
bool isSpecialNonInteracting(const T &p)
Identify a special non-interacting particles.
bool isFromHadron(T p, U hadron, bool &fromTau, bool &fromBSM)
Function to classify the particle.
bool isStable(const T &p)
Identify if the particle is stable, i.e. has not decayed.
bool isGeantino(const T &p)
bool isSimInteracting(const T &p)
Identify if the particle could interact with the detector during the simulation, e....
bool isEMInteracting(const T &p)
void findParticleAncestors(T thePart, std::set< T > &allancestors)
Function to find all ancestors of the particle.
bool isStrongInteracting(const T &p)
bool isInteracting(const T &p)
Identify if the particle with given PDG ID would not interact with the detector, i....
bool isMuon(const T &p)
bool isSUSY(const T &p)
auto findFinalStateParticles(V theVert) -> decltype(theVert->particles_out())
Function to find the stable particle descendants of the given vertex..
static const int NU_MU
bool isChargedNonShowering(const T &p)
Identify if the particle with given PDG ID would produce ID tracks but not shower in the detector if ...
bool isDecayed(const T &p)
Identify if the particle decayed.
T findMother(T thePart)
Function to get a mother of particle. MCTruthClassifier legacy.
auto particles_in(const HepMC::GenVertex *p)
void findParticleStableDescendants(T thePart, std::set< T > &allstabledescendants)
Function to get the particle stable MC daughters.
static const int NU_E
bool isBeam(const T &p)
Identify if the particle is beam particle.
static const int NU_TAU
bool isFinalState(const T &p)
Identify if the particle is final state particle.
bool isGenStable(const T &p)
Determine if the particle is stable at the generator (not det-sim) level,.
bool isHadron(const T &p)
bool isSimStable(const T &p)
Identify if the particle is considered stable at the post-detector-sim stage.
auto findSimulatedEndVertex(T thePart) -> decltype(thePart->end_vertex())
Function to find the end vertex of a particle.
bool isTau(const T &p)
bool isPhysical(const T &p)
Identify if the particle is physical, i.e. is stable or decayed.
bool isBSM(const T &p)
APID: graviton and all Higgs extensions are BSM.
Author: James Monk (jmonk@cern.ch)