ATLAS Offline Software
Functions
MC Namespace Reference

Functions

template<class T >
bool isInteracting (const T &p)
 Identify if the particle with given PDG ID would not interact with the detector, i.e. not a neutrino or WIMP. More...
 
template<class T >
bool isChargedNonShowering (const T &p)
 Identify if the particle with given PDG ID would produce ID tracks but not shower in the detector if stable. More...
 
template<class T >
bool isBeam (const T &p)
 
template<class T >
bool isDecayed (const T &p)
 
template<class T >
bool isStable (const T &p)
 
template<class T >
bool isFinalState (const T &p)
 
template<class T >
bool isPhysical (const T &p)
 
template<class T >
bool isPhysicalHadron (const T &p)
 
template<class T >
bool isGenStable (const T &p)
 Determine if the particle is stable at the generator (not det-sim) level,. More...
 
template<class T >
bool isSimStable (const T &p)
 Identify if the particle is considered stable at the post-detector-sim stage. More...
 
template<class T >
bool isSimInteracting (const T &p)
 Identify if the particle could interact with the detector during the simulation, e.g. not a neutrino or WIMP. More...
 
template<class T >
bool isStableOrSimDecayed (const T &p)
 Identify if particle is satble or decayed in simulation. More...
 
template<class T >
bool isZeroEnergyPhoton (const T &p)
 Identify a photon with zero energy. Probably a workaround for a generator bug. More...
 
template<class T >
bool isSingleParticle (const T &p)
 
template<class T >
bool isSpecialNonInteracting (const T &p)
 
template<class T >
getMother (T thePart)
 Function to get a mother of particle. MCTruthClassifier legacy. More...
 
template<class C , class T >
find_matching (C TruthTES, T bcin)
 Function to find a particle in container. More...
 
template<class T >
void findAllJetMothers (T thePart, std::set< T > &allJetMothers)
 Function to find all ancestors of the particle. More...
 
template<class T >
void findParticleDaughters (T thePart, std::set< T > &daughters)
 Function to get the particle stable MC daughters. More...
 
template<class T >
isHadronFromB (T p)
 Function to get the parent B hadron. More...
 
template<class T >
bool isHardScatVrtx (T pVert)
 Function to classify the vertex as hard scattering vertex. More...
 
template<class T >
bool fromHadron (T p, T hadptr, bool &fromTau, bool &fromBSM)
 Function to classify the particle. More...
 
template<class T >
auto findEndVert (T thePart) -> decltype(thePart->end_vertex())
 Function to find the end vertex of a particle. More...
 
template<class V >
auto findFinalStatePart (V EndVert) -> decltype(EndVert->particles_out())
 Function to find the stable particle descendants of the gived vertex.. More...
 

Function Documentation

◆ find_matching()

template<class C , class T >
T MC::find_matching ( TruthTES,
bcin 
)

Function to find a particle in container.

This can be used for HepMC3::GenVertexPtr, HepMC3::ConstGenVertexPtr or xAOD::TruthVertex*

Definition at line 103 of file HepMCHelpers.h.

103  {
104  T ptrPart = nullptr;
105  if (!bcin) return ptrPart;
106  for (T truthParticle : *TruthTES) {
107  if (HepMC::is_sim_descendant(bcin,truthParticle)) {
108  ptrPart = truthParticle;
109  break;
110  }
111  }
112  return ptrPart;
113  }

◆ findAllJetMothers()

template<class T >
void MC::findAllJetMothers ( thePart,
std::set< T > &  allJetMothers 
)

Function to find all ancestors of the particle.

This can be used for HepMC3::GenParticlePtr, HepMC3::ConstGenParticlePtr or xAOD::TruthParticle*

Definition at line 116 of file HepMCHelpers.h.

116  {
117  auto partOriVert = thePart->production_vertex();
118  if (!partOriVert) return;
119  auto incoming = partOriVert->particles_in();
120  for (auto theMoth: incoming) {
121  if (!theMoth) continue;
122  allJetMothers.insert(theMoth);
123  findAllJetMothers(theMoth, allJetMothers);
124  }
125  }

◆ findEndVert()

template<class T >
auto MC::findEndVert ( thePart) -> decltype(thePart->end_vertex())

Function to find the end vertex of a particle.

This algorithm allows for 1->1 decays.
This can be used for HepMC3::GenVertexPtr, HepMC3::ConstGenVertexPtr or xAOD::TruthVertex* and particle counterparts

Definition at line 212 of file HepMCHelpers.h.

212  {
213  decltype(thePart->end_vertex()) EndVert = thePart->end_vertex();
214  decltype(thePart->end_vertex()) pVert(nullptr);
215  if (EndVert != nullptr) {
216  do {
217  bool samePart = false;
218  pVert = nullptr;
219  auto outgoing = EndVert->particles_out();
220  auto incoming = EndVert->particles_in();
221  for (const auto& itrDaug: outgoing) {
222  if (!itrDaug) continue;
223  if (((itrDaug && HepMC::is_same_generator_particle(itrDaug,thePart)) ||
224  // brem on generator level for tau
225  (outgoing.size() == 1 && incoming.size() == 1 &&
227  itrDaug->pdg_id() == thePart->pdg_id()) {
228  samePart = true;
229  pVert = itrDaug->end_vertex();
230  }
231  }
232  if (samePart) EndVert = pVert;
233  } while (pVert != nullptr && pVert != EndVert); // pVert!=EndVert to prevent Sherpa loop
234  }
235  return EndVert;
236  }

◆ findFinalStatePart()

template<class V >
auto MC::findFinalStatePart ( EndVert) -> decltype(EndVert->particles_out())

Function to find the stable particle descendants of the gived vertex..

This can be used for HepMC3::GenVertexPtr, HepMC3::ConstGenVertexPtr or xAOD::TruthVertex* and particle counterparts

Definition at line 240 of file HepMCHelpers.h.

240  {
241  if (!EndVert) return {};
242  decltype(EndVert->particles_out()) finalStatePart;
243  auto outgoing = EndVert->particles_out();
244  for (const auto& thePart: outgoing) {
245  if (!thePart) continue;
246  finalStatePart.push_back(thePart);
247  if (isStable(thePart)) continue;
248  V pVert = findEndVert(thePart);
249  if (pVert == EndVert) break; // to prevent Sherpa loop
250  if (pVert != nullptr) {
251  auto vecPart = findFinalStatePart<V>(pVert);
252  finalStatePart.insert(finalStatePart.end(),vecPart.begin(),vecPart.end());
253  }
254  }
255  return finalStatePart;
256  }

◆ findParticleDaughters()

template<class T >
void MC::findParticleDaughters ( thePart,
std::set< T > &  daughters 
)

Function to get the particle stable MC daughters.

This can be used for HepMC3::GenParticlePtr, HepMC3::ConstGenParticlePtr or xAOD::TruthParticle*

Definition at line 129 of file HepMCHelpers.h.

129  {
130  auto endVtx = thePart->end_vertex();
131  if (!endVtx) return;
132  for (auto theDaughter: endVtx->particles_out()) {
133  if (theDaughter && isStable(theDaughter) && !HepMC::is_simulation_particle(theDaughter)) {
134  daughters.insert(theDaughter);
135  }
136  findParticleDaughters(theDaughter, daughters);
137  }
138  }

◆ fromHadron()

template<class T >
bool MC::fromHadron ( p,
hadptr,
bool &  fromTau,
bool &  fromBSM 
)

Function to classify the particle.

AV: This is MCtruthClassifier legacy. The function should be improved in the future. This can be used for HepMC3::GenVertexPtr, HepMC3::ConstGenVertexPtr or xAOD::TruthVertex*

Definition at line 184 of file HepMCHelpers.h.

184  {
185  if (isHadron(p)&&!isBeam(p)) return true; // trivial case
186  auto vtx = p->production_vertex();
187  if (!vtx) return false;
188  bool fromHad = false;
189  auto incoming = vtx->particles_in();
190  for (auto parent: incoming) {
191  if (!parent) continue;
192  // should this really go into parton-level territory?
193  // probably depends where BSM particles are being decayed
194  fromBSM |= isBSM(parent);
195  // sometimes Athena replaces status 2 with HepMC::SPECIALSTATUS, see e.g.
196  // PhysicsAnalysis/TruthParticleID/McParticleTools/src/EtaPtFilterTool.cxx#L374
197  // not at all clear why and unfortunately there's no documentation in the code
198  if (!isPhysical(parent) && parent->status() != HepMC::SPECIALSTATUS) return false;
199  fromTau |= isTau(parent);
200  if (isHadron(parent)&&!isBeam(parent)) {
201  if (!hadptr) hadptr = parent; // assumes linear hadron parentage
202  return true;
203  }
204  fromHad |= fromHadron(parent, hadptr, fromTau, fromBSM);
205  }
206  return fromHad;
207  }

◆ getMother()

template<class T >
T MC::getMother ( thePart)

Function to get a mother of particle. MCTruthClassifier legacy.

This can be used for HepMC3::GenVertexPtr, HepMC3::ConstGenVertexPtr or xAOD::TruthVertex*

Definition at line 68 of file HepMCHelpers.h.

68  {
69  auto partOriVert = thePart->production_vertex();
70  if (!partOriVert) return nullptr;
71 
72  long partPDG = thePart->pdg_id();
73  long MotherPDG(0);
74 
75  auto MothOriVert = partOriVert;
76  MothOriVert = nullptr;
77  T theMoth(nullptr);
78 
79  size_t itr = 0;
80  do {
81  if (itr != 0) partOriVert = MothOriVert;
82  auto incoming = partOriVert->particles_in();
83  for ( auto p: incoming) {
84  theMoth = p;
85  if (!theMoth) continue;
86  MotherPDG = theMoth->pdg_id();
87  MothOriVert = theMoth->production_vertex();
88  if (MotherPDG == partPDG) break;
89  }
90  itr++;
91  if (itr > 100) {
92  break;
93  }
94  } while (MothOriVert != nullptr && MotherPDG == partPDG && !HepMC::is_simulation_particle(thePart) &&
95  MothOriVert != partOriVert);
96  return theMoth;
97  }

◆ isBeam()

template<class T >
bool MC::isBeam ( const T &  p)
inline

Definition at line 28 of file HepMCHelpers.h.

28 { return p->status()%HepMC::SIM_STATUS_THRESHOLD == 4;}

◆ isChargedNonShowering()

template<class T >
bool MC::isChargedNonShowering ( const T &  p)
inline

Identify if the particle with given PDG ID would produce ID tracks but not shower in the detector if stable.

Definition at line 26 of file HepMCHelpers.h.

26 { return (isMuon<T>(p) || isSUSY<T>(p)); }

◆ isDecayed()

template<class T >
bool MC::isDecayed ( const T &  p)
inline

Definition at line 29 of file HepMCHelpers.h.

29 { return p->status()%HepMC::SIM_STATUS_THRESHOLD == 2;}

◆ isFinalState()

template<class T >
bool MC::isFinalState ( const T &  p)
inline

Definition at line 31 of file HepMCHelpers.h.

31 { return p->status()%HepMC::SIM_STATUS_THRESHOLD == 1 && !p->end_vertex();}

◆ isGenStable()

template<class T >
bool MC::isGenStable ( const T &  p)
inline

Determine if the particle is stable at the generator (not det-sim) level,.

Definition at line 36 of file HepMCHelpers.h.

36 { return isStable<T>(p) && !HepMC::is_simulation_particle<T>(p);}

◆ isHadronFromB()

template<class T >
T MC::isHadronFromB ( p)

Function to get the parent B hadron.

This can be used for HepMC3::GenParticlePtr, HepMC3::ConstGenParticlePtr or xAOD::TruthParticle*

AV: Strictly speaking, this is wrong and one has to check all the incoming particles. However that is MCTruthCalssifier legacy that should be fixed later at some point.

Definition at line 142 of file HepMCHelpers.h.

142  {
143  if (!p) return nullptr;
144  int pid = abs(p->pdg_id());
145  if (isBottomHadron(pid) || pid == BQUARK) return p;
146  if (pid == CQUARK || isNucleus(pid) || isBSM(pid) || !p->production_vertex()) return nullptr;
147  auto incoming = p->production_vertex()->particles_in();
148  if (incoming.size() == 0) return nullptr;
151  return isHadronFromB(incoming.front());
152  }

◆ isHardScatVrtx()

template<class T >
bool MC::isHardScatVrtx ( pVert)

Function to classify the vertex as hard scattering vertex.

AV: This is MCtruthClassifier legacy. Note that this function willnot capture some cases of the HardScattering vertices. The function should be improved in the future. This can be used for HepMC3::GenVertexPtr, HepMC3::ConstGenVertexPtr or xAOD::TruthVertex*

Definition at line 158 of file HepMCHelpers.h.

158  {
159  if (pVert == nullptr) return false;
160  T pV = pVert;
161  int numOfPartIn(0);
162  int pdg(0);
163 
164  do {
165  pVert = pV;
166  auto incoming = pVert->particles_in();
167  numOfPartIn = incoming.size();
168  pdg = numOfPartIn && incoming.front() != nullptr ? incoming.front()->pdg_id() : 0;
169  pV = numOfPartIn && incoming.front() != nullptr ? incoming.front()->production_vertex() : nullptr;
170 
171  } while (numOfPartIn == 1 && (std::abs(pdg) < 81 || std::abs(pdg) > 100) && pV != nullptr);
172 
173  if (numOfPartIn == 2) {
174  auto incoming = pVert->particles_in();
175  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;
176  }
177  return false;
178 }

◆ isInteracting()

template<class T >
bool MC::isInteracting ( const T &  p)
inline

Identify if the particle with given PDG ID would not interact with the detector, i.e. not a neutrino or WIMP.

Definition at line 23 of file HepMCHelpers.h.

23 { return isStrongInteracting<T>(p) || isEMInteracting<T>(p) || isGeantino<T>(p); }

◆ isPhysical()

template<class T >
bool MC::isPhysical ( const T &  p)
inline

Definition at line 32 of file HepMCHelpers.h.

32 { return isStable<T>(p) || isDecayed<T>(p); }

◆ isPhysicalHadron()

template<class T >
bool MC::isPhysicalHadron ( const T &  p)
inline

Definition at line 33 of file HepMCHelpers.h.

33 { return isHadron<T>(p) && isPhysical<T>(p);}

◆ isSimInteracting()

template<class T >
bool MC::isSimInteracting ( const T &  p)
inline

Identify if the particle could interact with the detector during the simulation, e.g. not a neutrino or WIMP.

Definition at line 42 of file HepMCHelpers.h.

42 { return isGenStable<T>(p) && isInteracting<T>(p);}

◆ isSimStable()

template<class T >
bool MC::isSimStable ( const T &  p)
inline

Identify if the particle is considered stable at the post-detector-sim stage.

Definition at line 39 of file HepMCHelpers.h.

39 { return isStable<T>(p) && !p->end_vertex() && HepMC::is_simulation_particle<T>(p);}

◆ isSingleParticle()

template<class T >
bool MC::isSingleParticle ( const T &  p)
inline

Definition at line 55 of file HepMCHelpers.h.

55 { return p->barcode() == HepMC::SINGLE_PARTICLE;}

◆ isSpecialNonInteracting()

template<class T >
bool MC::isSpecialNonInteracting ( const T &  p)
inline

Definition at line 57 of file HepMCHelpers.h.

57  {
58  const int apid = std::abs(p->pdg_id());
59  if (apid == 12 || apid == 14 || apid == 16) return true; //< neutrinos
60  if (apid == 1000022 || apid == 1000024 || apid == 5100022) return true; // SUSY & KK photon and Z partners
61  if (apid == 39 || apid == 1000039 || apid == 5000039) return true; //< gravitons: standard, SUSY and KK
62  if (apid == 9000001 || apid == 9000002 || apid == 9000003 || apid == 9000004 || apid == 9000005 || apid == 9000006) return true; //< exotic particles from monotop model
63  return false;
64  }

◆ isStable()

template<class T >
bool MC::isStable ( const T &  p)
inline

Definition at line 30 of file HepMCHelpers.h.

30 { return p->status()%HepMC::SIM_STATUS_THRESHOLD == 1;}

◆ isStableOrSimDecayed()

template<class T >
bool MC::isStableOrSimDecayed ( const T &  p)
inline

Identify if particle is satble or decayed in simulation.

  • a pathological case of decayed particle w/o end vertex. The decayed particles w/o end vertex might occur in case of simulation of long lived particles in Geant stripped off the decay products. I.e. those particles should be re-decayed later.

Definition at line 47 of file HepMCHelpers.h.

47  {
48  const auto vertex = p->end_vertex();
49  return ( isStable<T>(p) || (isDecayed<T>(p) && (!vertex || HepMC::is_simulation_vertex(vertex))));
50  }

◆ isZeroEnergyPhoton()

template<class T >
bool MC::isZeroEnergyPhoton ( const T &  p)
inline

Identify a photon with zero energy. Probably a workaround for a generator bug.

Definition at line 53 of file HepMCHelpers.h.

53 { return isPhoton<T>(p) && p->e() == 0;}
isNucleus
bool isNucleus(const T &p)
Definition: AtlasPID.h:212
python.PerfMonSerializer.p
def p
Definition: PerfMonSerializer.py:743
MC::fromHadron
bool fromHadron(T p, T hadptr, bool &fromTau, bool &fromBSM)
Function to classify the particle.
Definition: HepMCHelpers.h:184
isBSM
bool isBSM(const T &p)
Definition: AtlasPID.h:213
MC::findAllJetMothers
void findAllJetMothers(T thePart, std::set< T > &allJetMothers)
Function to find all ancestors of the particle.
Definition: HepMCHelpers.h:116
const
bool const RAWDATA *ch2 const
Definition: LArRodBlockPhysicsV0.cxx:562
isBottomHadron
bool isBottomHadron(const T &p)
Definition: AtlasPID.h:466
MC::isPhysical
bool isPhysical(const T &p)
Definition: HepMCHelpers.h:32
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:139
HepMC::SPECIALSTATUS
constexpr int SPECIALSTATUS
Constant that the meaning of which is currently lost, to be recovered...
Definition: MagicNumbers.h:39
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:129
MC::findEndVert
auto findEndVert(T thePart) -> decltype(thePart->end_vertex())
Function to find the end vertex of a particle.
Definition: HepMCHelpers.h:212
MC::isHadronFromB
T isHadronFromB(T p)
Function to get the parent B hadron.
Definition: HepMCHelpers.h:142
ParticleGun_EoverP_Config.pid
pid
Definition: ParticleGun_EoverP_Config.py:62
test_pyathena.parent
parent
Definition: test_pyathena.py:15
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:135
isTau
bool isTau(const T &p)
Definition: AtlasPID.h:148
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:36
HepMC::SINGLE_PARTICLE
constexpr int SINGLE_PARTICLE
Definition: MagicNumbers.h:53
isHadron
bool isHadron(const T &p)
Definition: AtlasPID.h:207
Trk::vertex
@ vertex
Definition: MeasurementType.h:21
MC::isStable
bool isStable(const T &p)
Definition: HepMCHelpers.h:30
MC::isBeam
bool isBeam(const T &p)
Definition: HepMCHelpers.h:28
if
if(febId1==febId2)
Definition: LArRodBlockPhysicsV0.cxx:569
MC::findParticleDaughters
void findParticleDaughters(T thePart, std::set< T > &daughters)
Function to get the particle stable MC daughters.
Definition: HepMCHelpers.h:129
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:143
TSU::T
unsigned long long T
Definition: L1TopoDataTypes.h:35