ATLAS Offline Software
HepMCHelpers.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2024 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>
13 
17 
18 namespace MC
19 {
20  namespace Pythia8
21  {
23  template <class T> inline bool isConditionA(const T& p) { return p->status() == 62 || p->status() == 52 || p->status() == 21 || p->status() == 22;}
24 
25  template <class T> inline bool isConditionB(const T& p) { return p->status() == 23;}
26 
27  template <class T> inline bool isConditionC(const T& p) { return p->status() > 30 && p->status() < 40;}
28  }
29 
30 #include "AtlasPID.h"
31 
33  template <class T> inline bool isInteracting(const T& p) { return isStrongInteracting<T>(p) || isEMInteracting<T>(p) || isGeantino<T>(p); }
34 
36  template <class T> inline bool isChargedNonShowering(const T& p) { return (isMuon<T>(p) || isSUSY<T>(p)); }
37 
39  template <class T> inline bool isBeam(const T& p) { return HepMC::status(p)%HepMC::SIM_STATUS_THRESHOLD == 4;}
40 
42  template <class T> inline bool isDecayed(const T& p) { return HepMC::status(p)%HepMC::SIM_STATUS_THRESHOLD == 2;}
43 
45  template <class T> inline bool isStable(const T& p) { return HepMC::status(p)%HepMC::SIM_STATUS_THRESHOLD == 1;}
46 
48  template <class T> inline bool isFinalState(const T& p) { return HepMC::status(p)%HepMC::SIM_STATUS_THRESHOLD == 1 && !p->end_vertex();}
49 
51  template <class T> inline bool isPhysical(const T& p) { return isStable<T>(p) || isDecayed<T>(p); }
52 
54  template <class T> inline bool isGenStable(const T& p) { return isStable<T>(p) && !HepMC::is_simulation_particle<T>(p);}
55 
57  template <class T> inline bool isSimStable(const T& p) { return isStable<T>(p) && !p->end_vertex() && HepMC::is_simulation_particle<T>(p);}
58 
60  template <class T> inline bool isSimInteracting(const T& p) { return isGenStable<T>(p) && isInteracting<T>(p);}
61 
65  template <class T> inline bool isStableOrSimDecayed(const T& p) {
66  const auto vertex = p->end_vertex();
67  return ( isStable<T>(p) || (isDecayed<T>(p) && (!vertex || HepMC::is_simulation_vertex(vertex))));
68  }
69 
71  template <class T> inline bool isZeroEnergyPhoton(const T& p) { return isPhoton<T>(p) && p->e() == 0;}
72 
74  template <class T> inline bool isSpecialNonInteracting(const T& p) {
75  const int apid = std::abs(p->pdg_id());
76  if (apid == NU_E || apid == NU_MU || apid == NU_TAU) return true; //< neutrinos
77  if (apid == 1000022 || apid == 1000024 || apid == 5100022) return true; // SUSY & KK photon and Z partners
78  if (apid == GRAVITON || apid == 1000039 || apid == 5000039) return true; //< gravitons: standard, SUSY and KK
79  if (apid == 9000001 || apid == 9000002 || apid == 9000003 || apid == 9000004 || apid == 9000005 || apid == 9000006) return true; //< exotic particles from monotop model
80  return false;
81  }
82 
84 
85  template <class T> T findMother(T thePart) {
86  auto partOriVert = thePart->production_vertex();
87  if (!partOriVert) return nullptr;
88 
89  long partPDG = thePart->pdg_id();
90  long MotherPDG(0);
91 
92  auto MothOriVert = partOriVert;
93  MothOriVert = nullptr;
94  T theMoth(nullptr);
95 
96  size_t itr = 0;
97  do {
98  if (itr != 0) partOriVert = MothOriVert;
99  auto incoming = partOriVert->particles_in();
100  for ( auto p: incoming) {
101  theMoth = p;
102  if (!theMoth) continue;
103  MotherPDG = theMoth->pdg_id();
104  MothOriVert = theMoth->production_vertex();
105  if (MotherPDG == partPDG) break;
106  }
107  itr++;
108  if (itr > 100) {
109  break;
110  }
111  } while (MothOriVert != nullptr && MotherPDG == partPDG && !HepMC::is_simulation_particle(thePart) &&
112  MothOriVert != partOriVert);
113  return theMoth;
114  }
115 
117 
118  template <class C, class T> T findMatching(C TruthContainer, T p) {
119  T ptrPart = nullptr;
120  if (!p) return ptrPart;
121  if constexpr (std::is_pointer_v<C> || HepMC::is_smart_ptr_v<C>){ //C is ptr
122  for (T truthParticle : *TruthContainer) {
123  if (HepMC::is_sim_descendant(p,truthParticle)) {
124  ptrPart = truthParticle;
125  break;
126  }
127  }
128  }
129  else {
130  for (T truthParticle : TruthContainer) {
131  if (HepMC::is_sim_descendant(p,truthParticle)) {
132  ptrPart = truthParticle;
133  break;
134  }
135  }
136  }
137  return ptrPart;
138  }
140 
141  template <class T> void findParticleAncestors(T thePart, std::set<T>& allancestors) {
142  auto prodVtx = thePart->production_vertex();
143  if (!prodVtx) return;
144  for (auto theMother: prodVtx->particles_in()) {
145  if (!theMother) continue;
146  allancestors.insert(theMother);
147  findParticleAncestors(theMother, allancestors);
148  }
149  }
150 
152 
153  template <class T> void findParticleStableDescendants(T thePart, std::set<T>& allstabledescendants) {
154  auto endVtx = thePart->end_vertex();
155  if (!endVtx) return;
156  for (auto theDaughter: endVtx->particles_out()) {
157  if (!theDaughter) continue;
158  if (isStable(theDaughter) && !HepMC::is_simulation_particle(theDaughter)) {
159  allstabledescendants.insert(theDaughter);
160  }
161  findParticleStableDescendants(theDaughter, allstabledescendants);
162  }
163  }
164 
168 
169  template <class T> bool isHardScatteringVertex(T pVert) {
170  if (pVert == nullptr) return false;
171  T pV = pVert;
172  int numOfPartIn(0);
173  int pdg(0);
174 
175  do {
176  pVert = pV;
177  auto incoming = pVert->particles_in();
178  numOfPartIn = incoming.size();
179  pdg = numOfPartIn && incoming.front() != nullptr ? incoming.front()->pdg_id() : 0;
180  pV = numOfPartIn && incoming.front() != nullptr ? incoming.front()->production_vertex() : nullptr;
181 
182  } while (numOfPartIn == 1 && (std::abs(pdg) < 81 || std::abs(pdg) > 100) && pV != nullptr);
183 
184  if (numOfPartIn == 2) {
185  auto incoming = pVert->particles_in();
186  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;
187  }
188  return false;
189 }
190 
194 
195  template <class T> bool isFromHadron(T p, T hadron, bool &fromTau, bool &fromBSM) {
196  if (isHadron(p)&&!isBeam(p)) return true; // trivial case
197  auto vtx = p->production_vertex();
198  if (!vtx) return false;
199  bool fromHad = false;
200  auto incoming = vtx->particles_in();
201  for (auto parent: incoming) {
202  if (!parent) continue;
203  // should this really go into parton-level territory?
204  // probably depends where BSM particles are being decayed
205  fromBSM |= isBSM(parent);
206  if (!isPhysical(parent)) return false;
207  fromTau |= isTau(parent);
208  if (isHadron(parent)&&!isBeam(parent)) {
209  if (!hadron) hadron = parent; // assumes linear hadron parentage
210  return true;
211  }
212  fromHad |= isFromHadron(parent, hadron, fromTau, fromBSM);
213  }
214  return fromHad;
215  }
216 
219 
220  template <class T> auto findSimulatedEndVertex(T thePart) -> decltype(thePart->end_vertex()) {
221  decltype(thePart->end_vertex()) EndVert = thePart->end_vertex();
222  decltype(thePart->end_vertex()) pVert(nullptr);
223  if (EndVert != nullptr) {
224  do {
225  bool samePart = false;
226  pVert = nullptr;
227  auto outgoing = EndVert->particles_out();
228  auto incoming = EndVert->particles_in();
229  for (const auto& itrDaug: outgoing) {
230  if (!itrDaug) continue;
231  if ((( HepMC::is_same_generator_particle(itrDaug,thePart)) ||
232  // brem on generator level for tau
233  (outgoing.size() == 1 && incoming.size() == 1 &&
235  itrDaug->pdg_id() == thePart->pdg_id()) {
236  samePart = true;
237  pVert = itrDaug->end_vertex();
238  }
239  }
240  if (samePart) EndVert = pVert;
241  } while (pVert != nullptr && pVert != EndVert); // pVert!=EndVert to prevent Sherpa loop
242  }
243  return EndVert;
244  }
245 
247 
248  template <class V> auto findFinalStateParticles(V theVert) -> decltype(theVert->particles_out()) {
249  if (!theVert) return {};
250  decltype(theVert->particles_out()) finalStatePart;
251  auto outgoing = theVert->particles_out();
252  for (const auto& thePart: outgoing) {
253  if (!thePart) continue;
254  finalStatePart.push_back(thePart);
255  if (isStable(thePart)) continue;
256  V pVert = findSimulatedEndVertex(thePart);
257  if (pVert == theVert) break; // to prevent Sherpa loop
258  if (pVert != nullptr) {
259  auto vecPart = findFinalStateParticles<V>(pVert);
260  finalStatePart.insert(finalStatePart.end(),vecPart.begin(),vecPart.end());
261  }
262  }
263  return finalStatePart;
264  }
265 
266 }
267 #endif
MC::isFromHadron
bool isFromHadron(T p, T hadron, bool &fromTau, bool &fromBSM)
Function to classify the particle.
Definition: HepMCHelpers.h:195
MCTruthPartClassifier::hadron
@ hadron
Definition: TruthClassifiers.h:148
MC::findParticleAncestors
void findParticleAncestors(T thePart, std::set< T > &allancestors)
Function to find all ancestors of the particle.
Definition: HepMCHelpers.h:141
MC::isSpecialNonInteracting
bool isSpecialNonInteracting(const T &p)
Identify a special non-interacting particles.
Definition: HepMCHelpers.h:74
isBSM
bool isBSM(const T &p)
APID: graviton and all Higgs extensions are BSM.
Definition: AtlasPID.h:841
MC::Pythia8::isConditionC
bool isConditionC(const T &p)
Definition: HepMCHelpers.h:27
MC
Definition: HepMCHelpers.h:19
MC::isGenStable
bool isGenStable(const T &p)
Determine if the particle is stable at the generator (not det-sim) level,.
Definition: HepMCHelpers.h:54
Pythia8
Author: James Monk (jmonk@cern.ch)
Definition: IPythia8Custom.h:10
MC::isPhysical
bool isPhysical(const T &p)
Identify if the particle is physical, i.e. is stable or decayed.
Definition: HepMCHelpers.h:51
python.utils.AtlRunQueryDQUtils.p
p
Definition: AtlRunQueryDQUtils.py:209
HepMC::is_same_generator_particle
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.
Definition: MagicNumbers.h:363
MC::findFinalStateParticles
auto findFinalStateParticles(V theVert) -> decltype(theVert->particles_out())
Function to find the stable particle descendants of the given vertex..
Definition: HepMCHelpers.h:248
MC::isFinalState
bool isFinalState(const T &p)
Identify if the particle is final state particle.
Definition: HepMCHelpers.h:48
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
MC::findMatching
T findMatching(C TruthContainer, T p)
Function to find a particle in container.
Definition: HepMCHelpers.h:118
test_pyathena.parent
parent
Definition: test_pyathena.py:15
MC::findParticleStableDescendants
void findParticleStableDescendants(T thePart, std::set< T > &allstabledescendants)
Function to get the particle stable MC daughters.
Definition: HepMCHelpers.h:153
HepMC::is_simulation_vertex
bool is_simulation_vertex(const T &v)
Method to establish if the vertex was created during simulation (TODO migrate to be based on status).
Definition: MagicNumbers.h:360
MC::findMother
T findMother(T thePart)
Function to get a mother of particle. MCTruthClassifier legacy.
Definition: HepMCHelpers.h:85
isTau
bool isTau(const T &p)
Definition: AtlasPID.h:205
MC::Pythia8::isConditionA
bool isConditionA(const T &p)
To be understood.
Definition: HepMCHelpers.h:23
HepMC::SIM_STATUS_THRESHOLD
constexpr int SIM_STATUS_THRESHOLD
Constant definiting the status threshold for simulated particles, eg. can be used to separate generat...
Definition: MagicNumbers.h:46
MagicNumbers.h
MC::isChargedNonShowering
bool isChargedNonShowering(const T &p)
Identify if the particle with given PDG ID would produce ID tracks but not shower in the detector if ...
Definition: HepMCHelpers.h:36
isHadron
bool isHadron(const T &p)
Definition: AtlasPID.h:348
MC::isStableOrSimDecayed
bool isStableOrSimDecayed(const T &p)
Identify if particle is satble or decayed in simulation.
Definition: HepMCHelpers.h:65
Trk::vertex
@ vertex
Definition: MeasurementType.h:21
MC::isStable
bool isStable(const T &p)
Identify if the particle is stable, i.e. has not decayed.
Definition: HepMCHelpers.h:45
MC::isDecayed
bool isDecayed(const T &p)
Identify if the particle decayed.
Definition: HepMCHelpers.h:42
MC::isBeam
bool isBeam(const T &p)
Identify if the particle is beam particle.
Definition: HepMCHelpers.h:39
MC::isSimInteracting
bool isSimInteracting(const T &p)
Identify if the particle could interact with the detector during the simulation, e....
Definition: HepMCHelpers.h:60
MC::Pythia8::isConditionB
bool isConditionB(const T &p)
Definition: HepMCHelpers.h:25
MC::isInteracting
bool isInteracting(const T &p)
Identify if the particle with given PDG ID would not interact with the detector, i....
Definition: HepMCHelpers.h:33
HepMC::is_sim_descendant
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,...
Definition: MagicNumbers.h:372
MC::isSimStable
bool isSimStable(const T &p)
Identify if the particle is considered stable at the post-detector-sim stage.
Definition: HepMCHelpers.h:57
MC::isZeroEnergyPhoton
bool isZeroEnergyPhoton(const T &p)
Identify a photon with zero energy. Probably a workaround for a generator bug.
Definition: HepMCHelpers.h:71
HepMC::status
int status(const T &p)
Definition: MagicNumbers.h:143
AtlasPID.h
MC::isHardScatteringVertex
bool isHardScatteringVertex(T pVert)
Function to classify the vertex as hard scattering vertex.
Definition: HepMCHelpers.h:169
MC::findSimulatedEndVertex
auto findSimulatedEndVertex(T thePart) -> decltype(thePart->end_vertex())
Function to find the end vertex of a particle.
Definition: HepMCHelpers.h:220