|
ATLAS Offline Software
|
Go to the documentation of this file.
10 namespace TruthDecoratorHelpers {
13 return particle_A->
pt() < particle_B->
pt();
18 if ( not truth ) {
return nullptr; }
22 if ( not truth_vertex || truth_vertex->
perp() > 440.0 ) {
30 if ( !vertex_A || !vertex_B ) {
return 999.0; }
31 return (vertex_A->
v4().Vect() - vertex_B->
v4().Vect()).Mag();
35 if( truth_particle ==
nullptr ) {
return false; }
36 if( flavour == 5 && truth_particle->
isBottomHadron() ) {
return true; }
37 if( flavour == 4 && truth_particle->
isCharmHadron() ) {
return true; }
42 if (!
is_bc_hadron(truth_particle, flavour))
return false;
57 if ( truth_particle ==
nullptr ) {
return nullptr; }
59 if (
depth>30) {
return nullptr; }
62 return truth_particle;
64 for (
unsigned int p = 0;
p < truth_particle->
nParents();
p++) {
66 if(
parent == truth_particle)
continue;
68 if ( parent_hadron !=
nullptr ) {
85 return TruthType::Label::Other;
98 return TruthSource::Label::Other;
103 std::vector<const xAOD::TruthVertex*>& seen_vertices,
104 const float truthVertexMergeDistance) {
110 if (
get_distance(this_vertex, truth_PV) < truthVertexMergeDistance) {
114 for (
size_t i = 0;
i != seen_vertices.size();
i++) {
116 if (
dr < truthVertexMergeDistance ) {
121 seen_vertices.push_back(this_vertex);
122 return seen_vertices.size();
const TruthParticle_v1 * parent(size_t i=0) const
Retrieve the i-th mother (TruthParticle) of this TruthParticle.
bool is_bc_hadron(const xAOD::TruthParticle *truth_particle, int flavour)
@ depth
pointing depth of the shower as calculated in egammaqgcld
Electron_v1 Electron
Definition of the current "egamma version".
bool isElectron() const
Whether the particle is an electron (or positron)
float get_distance(const xAOD::TruthVertex *vertex_A, const xAOD::TruthVertex *vertex_B)
bool isHadronicInteraction(int origin)
from hadronic interactions
This file contains "getter" functions used for accessing tagger inputs from the EDM.
bool isCharmHadron() const
Determine if the PID is that of a c-hadron.
bool is_weakly_decaying_hadron(const xAOD::TruthParticle *truth_particle)
bool isStrangeMesonDecay(int origin)
from strange meson decay
const xAOD::TruthParticle * get_parent_hadron(const xAOD::TruthParticle *truth_particle, bool user_called=true, int depth=0)
int get_source_type(const int origin)
Class providing the definition of the 4-vector interface.
FourVec_t v4() const
The full 4-vector of the vertex.
bool isGammaConversion(int origin)
from conversions
bool hasDecayVtx() const
Check for a decay vertex on this particle.
size_t nParents() const
Number of parents of this particle.
Class describing a truth particle in the MC record.
const xAOD::TruthVertex * get_truth_vertex(const xAOD::TruthParticle *truth)
float perp() const
Vertex transverse distance from the beam line.
virtual double pt() const =0
The transverse momentum ( ) of the particle.
const TruthVertex_v1 * decayVtx() const
The decay vertex of this particle.
bool isBottomHadron() const
Determine if the PID is that of a b-hadron.
int get_vertex_index(const xAOD::TruthVertex *vertex, const xAOD::TruthVertex *truth_PV, std::vector< const xAOD::TruthVertex * > &seen_vertices, const float truthVertexMergeDistance)
Class describing a truth vertex in the MC record.
bool isSecondary(int origin)
from long living particle decays or gamma conversions or hadronic interactions and anything else with...
std::vector< const TruthParticle * > particles_out() const
Get the outgoing particles.
Photon_v1 Photon
Definition of the current "egamma version".
struct TBPatternUnitContext Muon
bool isStrangeBaryonDecay(int origin)
from strange baryon decay
bool isPhoton() const
Whether the particle is a photon.
bool sort_particles(const xAOD::IParticle *particle_A, const xAOD::IParticle *particle_B)
bool isMuon() const
Whether the particle is a muon (or antimuon)
int get_truth_type(const xAOD::TruthParticle *truth_particle)
int pdgId() const
PDG ID code.
bool isStrangeMeson() const
Determine if the PID is that of a strange meson.
double charge() const
Physical charge.