ATLAS Offline Software
Loading...
Searching...
No Matches
DerivationFramework::TruthDecayCollectionMaker Class Reference

#include <TruthDecayCollectionMaker.h>

Inheritance diagram for DerivationFramework::TruthDecayCollectionMaker:
Collaboration diagram for DerivationFramework::TruthDecayCollectionMaker:

Public Member Functions

virtual StatusCode initialize () override final
virtual StatusCode addBranches (const EventContext &ctx) const override final

Private Member Functions

int addTruthParticle (const EventContext &ctx, const xAOD::TruthParticle &old_part, xAOD::TruthParticleContainer *part_cont, xAOD::TruthVertexContainer *vert_cont, std::vector< int > &seen_particles, const int generations=-1) const
int addTruthVertex (const EventContext &, const xAOD::TruthVertex &old_vert, xAOD::TruthParticleContainer *part_cont, xAOD::TruthVertexContainer *vert_cont, std::vector< int > &seen_particles, const int generations=-1) const
bool id_ok (const xAOD::TruthParticle &part) const

Private Attributes

Gaudi::Property< std::vector< int > > m_pdgIdsToKeep {this, "PDGIDsToKeep", {}, "PDG IDs of particles to build the collection from"}
 < List of PDG IDs to build this collection from
Gaudi::Property< bool > m_keepBHadrons {this, "KeepBHadrons", false, "Keep b-hadrons (easier than by PDG ID)"}
 < Option to keep all b-hadrons (better than giving PDG IDs)
Gaudi::Property< bool > m_keepCHadrons {this, "KeepCHadrons", false, "Keep c-hadrons (easier than by PDG ID)"}
 < Option to keep all c-hadrons (better than giving PDG IDs)
Gaudi::Property< bool > m_keepBSM {this, "KeepBSM", false, "Keep BSM particles (easier than by PDG ID)"}
 < Option to keep all BSM particles (better than giving PDG IDs)
Gaudi::Property< bool > m_rejectHadronChildren {this, "RejectHadronChildren", false, "Drop hadron descendants"}
 < Option to reject hadron descendants
SG::ReadHandleKey< xAOD::TruthParticleContainerm_particlesKey {this, "ParticlesKey", "TruthParticles", "ReadHandleKey for input TruthParticleContainer"}
SG::ReadDecorHandleKey< xAOD::TruthParticleContainerm_originAccessorKey {this, "Input_classifierParticleOrigin", m_particlesKey, "classifierParticleOrigin","Name of the decoration which records the particle origin as determined by the MCTruthClassifier"}
SG::ReadDecorHandleKey< xAOD::TruthParticleContainerm_typeAccessorKey {this, "Input_classifierParticleType", m_particlesKey, "classifierParticleType","Name of the decoration which records the particle type as determined by the MCTruthClassifier"}
SG::ReadDecorHandleKey< xAOD::TruthParticleContainerm_outcomeAccessorKey {this, "Input_classifierParticleOutCome", m_particlesKey, "classifierParticleOutCome","Name of the decoration which records the particle outcome as determined by the MCTruthClassifier"}
SG::ReadDecorHandleKey< xAOD::TruthParticleContainerm_classificationAccessorKey {this, "Input_Classification", m_particlesKey, "Classification","Name of the decoration which records the particle outcome as determined by the MCTruthClassifier"}
SG::WriteHandleKey< xAOD::TruthVertexContainerm_outputVerticesKey {this, "NewVertexKey", "", "WriteHandleKey for new TruthVertexContainer"}
SG::WriteHandleKey< xAOD::TruthParticleContainerm_outputParticlesKey {this, "NewParticleKey", "", "WriteHandleKey for new TruthParticleContainer"}
SG::WriteDecorHandleKey< xAOD::TruthParticleContainerm_originDecoratorKey {this, "classifierParticleOrigin", m_outputParticlesKey, "classifierParticleOrigin","Name of the decoration which records the particle origin as determined by the MCTruthClassifier"}
 FIXME Using WriteDecorHandles for decorations on a Container created in the current algorithm is unnecessary.
SG::WriteDecorHandleKey< xAOD::TruthParticleContainerm_typeDecoratorKey {this, "classifierParticleType", m_outputParticlesKey, "classifierParticleType","Name of the decoration which records the particle type as determined by the MCTruthClassifier"}
SG::WriteDecorHandleKey< xAOD::TruthParticleContainerm_outcomeDecoratorKey {this, "classifierParticleOutCome", m_outputParticlesKey, "classifierParticleOutCome","Name of the decoration which records the particle outcome as determined by the MCTruthClassifier"}
SG::WriteDecorHandleKey< xAOD::TruthParticleContainerm_classificationDecoratorKey {this, "Classification", m_outputParticlesKey, "Classification","Name of the decoration which records the particle outcome as determined by the MCTruthClassifier"}
SG::WriteDecorHandleKey< xAOD::TruthParticleContainerm_motherIDDecoratorKey {this, "motherID", m_outputParticlesKey, "motherID","Name of the decoration which records the ID of the particle's mother"}
SG::WriteDecorHandleKey< xAOD::TruthParticleContainerm_daughterIDDecoratorKey {this, "daughterID", m_outputParticlesKey, "daughterID","Name of the decoration which records the ID of the particle's daughter"}
Gaudi::Property< int > m_generations {this, "Generations", -1, "Number of generations after the particle in question to keep (-1 for all)"}
 ^^^^ These should be replaced by SG::Accessor.

Detailed Description

Definition at line 26 of file TruthDecayCollectionMaker.h.

Member Function Documentation

◆ addBranches()

StatusCode DerivationFramework::TruthDecayCollectionMaker::addBranches ( const EventContext & ctx) const
finaloverridevirtual

Definition at line 68 of file TruthDecayCollectionMaker.cxx.

69{
70 // Event context for AthenaMT
71
72 // Retrieve truth collections
73 SG::ReadHandle<xAOD::TruthParticleContainer> truthParticles(m_particlesKey,ctx);
74 if (!truthParticles.isValid()) {
75 ATH_MSG_ERROR("Couldn't retrieve TruthParticle collection with name " << m_particlesKey);
76 return StatusCode::FAILURE;
77 }
78
79 // Create the new particle containers and WriteHandles
80 SG::WriteHandle<xAOD::TruthParticleContainer> newParticlesWriteHandle(m_outputParticlesKey, ctx);
81 ATH_CHECK(newParticlesWriteHandle.record(std::make_unique<xAOD::TruthParticleContainer>(),
82 std::make_unique<xAOD::TruthParticleAuxContainer>()));
83 ATH_MSG_DEBUG( "Recorded new TruthParticleContainer with key: " << (m_outputParticlesKey.key()));
84 // Create the new vertex containers and WriteHandles
85 SG::WriteHandle<xAOD::TruthVertexContainer> newVerticesWriteHandle(m_outputVerticesKey, ctx);
86 ATH_CHECK(newVerticesWriteHandle.record(std::make_unique<xAOD::TruthVertexContainer>(),
87 std::make_unique<xAOD::TruthVertexAuxContainer>()));
88 ATH_MSG_DEBUG( "Recorded new TruthVertexContainer with key: " << (m_outputVerticesKey.key()));
89
90 // List of unique IDs for particles in our collection already. Because of the way we recurse,
91 // adding more particles as we go, there should be no need to add (or help from adding) the
92 // unique IDs of particles that we are *not* going to keep
93 std::vector<int> seen_particles;
94 // Go through that list of particles!
95 for (const auto * part : *truthParticles){
96 // If this passes my cuts, keep it
97 if (id_ok(*part)){
98 addTruthParticle( ctx, *part, newParticlesWriteHandle.ptr(), newVerticesWriteHandle.ptr(), seen_particles , m_generations );
99 }
100 } // Loop over the initial truth particle collection
101 return StatusCode::SUCCESS;
102}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_DEBUG(x)
SG::WriteHandleKey< xAOD::TruthVertexContainer > m_outputVerticesKey
Gaudi::Property< int > m_generations
^^^^ These should be replaced by SG::Accessor.
SG::WriteHandleKey< xAOD::TruthParticleContainer > m_outputParticlesKey
int addTruthParticle(const EventContext &ctx, const xAOD::TruthParticle &old_part, xAOD::TruthParticleContainer *part_cont, xAOD::TruthVertexContainer *vert_cont, std::vector< int > &seen_particles, const int generations=-1) const
SG::ReadHandleKey< xAOD::TruthParticleContainer > m_particlesKey
bool id_ok(const xAOD::TruthParticle &part) const

◆ addTruthParticle()

int DerivationFramework::TruthDecayCollectionMaker::addTruthParticle ( const EventContext & ctx,
const xAOD::TruthParticle & old_part,
xAOD::TruthParticleContainer * part_cont,
xAOD::TruthVertexContainer * vert_cont,
std::vector< int > & seen_particles,
const int generations = -1 ) const
private

Definition at line 104 of file TruthDecayCollectionMaker.cxx.

109 {
110 // See if we've seen it - note, could also do this with a unary function on the container itself
111 if (std::find(seen_particles.begin(),seen_particles.end(),HepMC::uniqueID(&old_part))!=seen_particles.end()){
112 for (size_t p=0;p<part_cont->size();++p){
113 // Was it a hit?
114 const xAOD::TruthParticle *theParticle = (*part_cont)[p];
115 if (HepMC::is_same_particle(theParticle,&old_part)) return p;
116 } // Look through the old container
117 } // Found it in the old container
118 // Now we have seen it
119 seen_particles.push_back(HepMC::uniqueID(&old_part));
120 // Set up decorators
121 SG::WriteDecorHandle<xAOD::TruthParticleContainer, unsigned int > originDecorator(m_originDecoratorKey, ctx);
122 SG::WriteDecorHandle<xAOD::TruthParticleContainer, unsigned int > typeDecorator(m_typeDecoratorKey, ctx);
123 SG::WriteDecorHandle<xAOD::TruthParticleContainer, unsigned int > outcomeDecorator(m_outcomeDecoratorKey, ctx);
124 SG::WriteDecorHandle<xAOD::TruthParticleContainer, unsigned int > classificationDecorator(m_classificationDecoratorKey, ctx);
125 SG::WriteDecorHandle< xAOD::TruthParticleContainer, int > motherIDDecorator(m_motherIDDecoratorKey, ctx);
126 SG::WriteDecorHandle< xAOD::TruthParticleContainer, int > daughterIDDecorator(m_daughterIDDecoratorKey, ctx);
127 // Make a truth particle and add it to the container
128 xAOD::TruthParticle* xTruthParticle = new xAOD::TruthParticle();
129 part_cont->push_back( xTruthParticle );
130 // Fill with numerical content
131 xTruthParticle->setPdgId(old_part.pdgId());
132 xTruthParticle->setUid(HepMC::uniqueID(&old_part));
133 xTruthParticle->setStatus(old_part.status());
134 xTruthParticle->setM(old_part.m());
135 xTruthParticle->setPx(old_part.px());
136 xTruthParticle->setPy(old_part.py());
137 xTruthParticle->setPz(old_part.pz());
138 xTruthParticle->setE(old_part.e());
139 // Copy over the polarization information if it's there
140 if (old_part.polarization().valid()){
143 }
144 // Make a link to this particle
145 int my_index = part_cont->size()-1;
146 ElementLink<xAOD::TruthParticleContainer> eltp(*part_cont, my_index);
147 // Decay vertex information
148 if (old_part.hasDecayVtx()) {
149 int vert_index = addTruthVertex( ctx, *old_part.decayVtx(), part_cont, vert_cont, seen_particles, generations);
150 ElementLink<xAOD::TruthVertexContainer> eltv( *vert_cont, vert_index );
151 xTruthParticle->setDecayVtxLink( eltv );
152 (*vert_cont)[vert_index]->addIncomingParticleLink( eltp );
153 }
154 // Copy over the decorations if they are available
155 SG::ReadDecorHandle<xAOD::TruthParticleContainer, unsigned int > classifierParticleTypeAcc(m_typeAccessorKey, ctx);
156 SG::ReadDecorHandle<xAOD::TruthParticleContainer, unsigned int > classifierParticleOriginAcc(m_originAccessorKey, ctx);
157 SG::ReadDecorHandle<xAOD::TruthParticleContainer, unsigned int > classifierParticleOutcomeAcc(m_outcomeAccessorKey, ctx);
158 SG::ReadDecorHandle<xAOD::TruthParticleContainer, unsigned int > ClassificationAcc(m_classificationAccessorKey, ctx);
159
160 typeDecorator(*xTruthParticle) = classifierParticleTypeAcc(old_part);
161 originDecorator(*xTruthParticle) = classifierParticleOriginAcc(old_part);
162 outcomeDecorator(*xTruthParticle) = classifierParticleOutcomeAcc(old_part);
163 classificationDecorator(*xTruthParticle) = ClassificationAcc(old_part);
164
165 // Return a link to this particle
166 return my_index;
167}
value_type push_back(value_type pElem)
Add an element to the end of the collection.
size_type size() const noexcept
Returns the number of elements in the collection.
SG::WriteDecorHandleKey< xAOD::TruthParticleContainer > m_daughterIDDecoratorKey
int addTruthVertex(const EventContext &, const xAOD::TruthVertex &old_vert, xAOD::TruthParticleContainer *part_cont, xAOD::TruthVertexContainer *vert_cont, std::vector< int > &seen_particles, const int generations=-1) const
SG::WriteDecorHandleKey< xAOD::TruthParticleContainer > m_outcomeDecoratorKey
SG::ReadDecorHandleKey< xAOD::TruthParticleContainer > m_outcomeAccessorKey
SG::WriteDecorHandleKey< xAOD::TruthParticleContainer > m_motherIDDecoratorKey
SG::ReadDecorHandleKey< xAOD::TruthParticleContainer > m_classificationAccessorKey
SG::ReadDecorHandleKey< xAOD::TruthParticleContainer > m_originAccessorKey
SG::ReadDecorHandleKey< xAOD::TruthParticleContainer > m_typeAccessorKey
SG::WriteDecorHandleKey< xAOD::TruthParticleContainer > m_typeDecoratorKey
SG::WriteDecorHandleKey< xAOD::TruthParticleContainer > m_originDecoratorKey
FIXME Using WriteDecorHandles for decorations on a Container created in the current algorithm is unne...
SG::WriteDecorHandleKey< xAOD::TruthParticleContainer > m_classificationDecoratorKey
void setPy(float value)
Set the y component of the particle's momentum.
int status() const
Status code.
void setUid(int value)
Set unique ID.
void setM(float value)
Also store the mass.
virtual double m() const override final
The mass of the particle.
void setE(float value)
Set the energy of the particle.
int pdgId() const
PDG ID code.
const TruthVertex_v1 * decayVtx() const
The decay vertex of this particle.
Polarization polarization() const
Retrieve a full Polarization with a single call.
void setPz(float value)
Set the z component of the particle's momentum.
float px() const
The x component of the particle's momentum.
bool hasDecayVtx() const
Check for a decay vertex on this particle.
virtual double e() const override final
The total energy of the particle.
bool setPolarizationParameter(float value, PolParam parameter)
Set method for polarization parameter values.
void setDecayVtxLink(const ElementLink< TruthVertexContainer > &link)
Set the decay vertex of the particle.
float py() const
The y component of the particle's momentum.
void setStatus(int value)
Set status code.
void setPdgId(int pid)
Set PDG ID code.
void setPx(float value)
Set the x component of the particle's momentum.
bool polarizationParameter(float &value, PolParam parameter) const
Accessor for polarization parameters.
float pz() const
The z component of the particle's momentum.
@ polarizationPhi
Polarization in ( )
@ polarizationTheta
Polarization in ( )
int uniqueID(const T &p)
bool is_same_particle(const T1 &p1, const T2 &p2)
Method to establish if two particles in the GenEvent actually represent the same particle.
TruthParticle_v1 TruthParticle
Typedef to implementation.
bool valid() const
Check if the stored values are valid.

◆ addTruthVertex()

int DerivationFramework::TruthDecayCollectionMaker::addTruthVertex ( const EventContext & ctx,
const xAOD::TruthVertex & old_vert,
xAOD::TruthParticleContainer * part_cont,
xAOD::TruthVertexContainer * vert_cont,
std::vector< int > & seen_particles,
const int generations = -1 ) const
private

Definition at line 169 of file TruthDecayCollectionMaker.cxx.

174 {
175 // Make a new vertex and add it to the container
176 xAOD::TruthVertex* xTruthVertex = new xAOD::TruthVertex();
177 vert_cont->push_back( xTruthVertex );
178 // Get a link to this vertex -- will be used to set production vertices on all the next particles
179 int my_index = vert_cont->size()-1;
180 ElementLink<xAOD::TruthVertexContainer> eltv(*vert_cont, my_index);
181 // Set properties
182 xTruthVertex->setStatus(HepMC::status(old_vert));
183 xTruthVertex->setUid(HepMC::uniqueID(old_vert));
184 xTruthVertex->setX(old_vert.x());
185 xTruthVertex->setY(old_vert.y());
186 xTruthVertex->setZ(old_vert.z());
187 xTruthVertex->setT(old_vert.t());
188 // If we are done, then stop here
189 if (generations==0) return my_index;
190 // Add all the outgoing particles
191 for (size_t n=0;n<old_vert.nOutgoingParticles();++n){
192 if (!old_vert.outgoingParticle(n)) continue; // Just in case we removed some truth particles, e.g. G4 decays
193 if (m_rejectHadronChildren && old_vert.outgoingParticle(n)->isHadron()) { // Option to skip hadrons outright{
194 continue;
195 }
196 // Continue on the next generation; note that we only decrement the generation if this particle doesn't also pass our cuts
197 int part_index = addTruthParticle( ctx, *old_vert.outgoingParticle(n), part_cont, vert_cont, seen_particles,
198 generations-1+(id_ok(*old_vert.outgoingParticle(n))?1:0) );
199 ElementLink<xAOD::TruthParticleContainer> eltp( *part_cont, part_index);
200 xTruthVertex->addOutgoingParticleLink( eltp );
201 (*part_cont)[part_index]->setProdVtxLink( eltv );
202 }
203 // Return a link to this vertex
204 return my_index;
205}
Gaudi::Property< bool > m_rejectHadronChildren
< Option to reject hadron descendants
bool isHadron() const
Whether the particle is a hadron.
const TruthParticle_v1 * outgoingParticle(size_t index) const
Get one of the outgoing particles.
void addOutgoingParticleLink(const TPLink_t &link)
Add one outgoing particle.
float z() const
Vertex longitudinal distance along the beam line form the origin.
void setStatus(int value)
Set the vertex status.
void setZ(float value)
Set the vertex's longitudinal distance from the origin.
float y() const
Vertex y displacement.
void setUid(int value)
Set the vertex unique ID.
float t() const
Vertex time.
size_t nOutgoingParticles() const
Get the number of outgoing particles.
void setT(float value)
Set the vertex time.
float x() const
Vertex x displacement.
void setX(float value)
Set the x displacement of the vertex.
void setY(float value)
Set the y displacement of the vertex.
int status(const T &p)
TruthVertex_v1 TruthVertex
Typedef to implementation.
Definition TruthVertex.h:15

◆ id_ok()

bool DerivationFramework::TruthDecayCollectionMaker::id_ok ( const xAOD::TruthParticle & part) const
private

Definition at line 207 of file TruthDecayCollectionMaker.cxx.

208{
209 // Check list of PDG IDs to keep
210 for (int id : m_pdgIdsToKeep){
211 if (part.absPdgId()==id){
212 return true;
213 } // Found a particle of interest!
214 } // Loop over the PDG IDs we want to keep
215 // Also check functions for B/C/BSM
216 return (m_keepBHadrons && part.isBottomHadron()) ||
217
218 (m_keepCHadrons && part.isCharmHadron()) ||
219
220 (m_keepBSM && part.isBSM());
221}
Gaudi::Property< bool > m_keepCHadrons
< Option to keep all c-hadrons (better than giving PDG IDs)
Gaudi::Property< std::vector< int > > m_pdgIdsToKeep
< List of PDG IDs to build this collection from
Gaudi::Property< bool > m_keepBHadrons
< Option to keep all b-hadrons (better than giving PDG IDs)
Gaudi::Property< bool > m_keepBSM
< Option to keep all BSM particles (better than giving PDG IDs)

◆ initialize()

StatusCode DerivationFramework::TruthDecayCollectionMaker::initialize ( )
finaloverridevirtual

Definition at line 30 of file TruthDecayCollectionMaker.cxx.

31{
32 ATH_MSG_VERBOSE("initialize() ...");
33
34 // Input truth particles
35 ATH_CHECK( m_particlesKey.initialize() );
36 ATH_MSG_INFO("Using " << m_particlesKey.key() << " as the input truth container key");
37
38 // Accessors (ReadDecorHandleKeys) TODO Convert these to SG::ConstAccessor?
39 ATH_CHECK(m_originAccessorKey.initialize());
40 ATH_CHECK(m_typeAccessorKey.initialize());
41 ATH_CHECK(m_outcomeAccessorKey.initialize());
43
44 // Output particle/vertex containers
45 ATH_CHECK(m_outputParticlesKey.initialize());
46 ATH_MSG_INFO("New truth particles container key: " << m_outputParticlesKey.key() );
47 ATH_CHECK(m_outputVerticesKey.initialize());
48 ATH_MSG_INFO("New truth vertices container key: " << m_outputVerticesKey.key() );
49
50 if (m_pdgIdsToKeep.empty() && !m_keepBHadrons && !m_keepCHadrons && !m_keepBSM) {
51 ATH_MSG_FATAL("No PDG IDs provided, not keeping b- or c-hadrons or BSM particles -- what do you want?");
52 return StatusCode::FAILURE;
53 }
54
55 // Decorators TODO Convert these to SG::Accessor?
56 ATH_CHECK(m_originDecoratorKey.initialize());
57 ATH_CHECK(m_typeDecoratorKey.initialize());
58 ATH_CHECK(m_outcomeDecoratorKey.initialize());
62
63 return StatusCode::SUCCESS;
64}
#define ATH_MSG_FATAL(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)

Member Data Documentation

◆ m_classificationAccessorKey

SG::ReadDecorHandleKey<xAOD::TruthParticleContainer> DerivationFramework::TruthDecayCollectionMaker::m_classificationAccessorKey {this, "Input_Classification", m_particlesKey, "Classification","Name of the decoration which records the particle outcome as determined by the MCTruthClassifier"}
private

Definition at line 55 of file TruthDecayCollectionMaker.h.

56{this, "Input_Classification", m_particlesKey, "Classification","Name of the decoration which records the particle outcome as determined by the MCTruthClassifier"};

◆ m_classificationDecoratorKey

SG::WriteDecorHandleKey<xAOD::TruthParticleContainer> DerivationFramework::TruthDecayCollectionMaker::m_classificationDecoratorKey {this, "Classification", m_outputParticlesKey, "Classification","Name of the decoration which records the particle outcome as determined by the MCTruthClassifier"}
private

Definition at line 70 of file TruthDecayCollectionMaker.h.

71{this, "Classification", m_outputParticlesKey, "Classification","Name of the decoration which records the particle outcome as determined by the MCTruthClassifier"};

◆ m_daughterIDDecoratorKey

SG::WriteDecorHandleKey<xAOD::TruthParticleContainer> DerivationFramework::TruthDecayCollectionMaker::m_daughterIDDecoratorKey {this, "daughterID", m_outputParticlesKey, "daughterID","Name of the decoration which records the ID of the particle's daughter"}
private

Definition at line 74 of file TruthDecayCollectionMaker.h.

75{this, "daughterID", m_outputParticlesKey, "daughterID","Name of the decoration which records the ID of the particle's daughter"};

◆ m_generations

Gaudi::Property<int> DerivationFramework::TruthDecayCollectionMaker::m_generations {this, "Generations", -1, "Number of generations after the particle in question to keep (-1 for all)"}
private

^^^^ These should be replaced by SG::Accessor.

< Number of generations after the particle in question to keep

Definition at line 78 of file TruthDecayCollectionMaker.h.

79{this, "Generations", -1, "Number of generations after the particle in question to keep (-1 for all)"};

◆ m_keepBHadrons

Gaudi::Property<bool> DerivationFramework::TruthDecayCollectionMaker::m_keepBHadrons {this, "KeepBHadrons", false, "Keep b-hadrons (easier than by PDG ID)"}
private

< Option to keep all b-hadrons (better than giving PDG IDs)

Definition at line 37 of file TruthDecayCollectionMaker.h.

38{this, "KeepBHadrons", false, "Keep b-hadrons (easier than by PDG ID)"};

◆ m_keepBSM

Gaudi::Property<bool> DerivationFramework::TruthDecayCollectionMaker::m_keepBSM {this, "KeepBSM", false, "Keep BSM particles (easier than by PDG ID)"}
private

< Option to keep all BSM particles (better than giving PDG IDs)

Definition at line 41 of file TruthDecayCollectionMaker.h.

42{this, "KeepBSM", false, "Keep BSM particles (easier than by PDG ID)"};

◆ m_keepCHadrons

Gaudi::Property<bool> DerivationFramework::TruthDecayCollectionMaker::m_keepCHadrons {this, "KeepCHadrons", false, "Keep c-hadrons (easier than by PDG ID)"}
private

< Option to keep all c-hadrons (better than giving PDG IDs)

Definition at line 39 of file TruthDecayCollectionMaker.h.

40{this, "KeepCHadrons", false, "Keep c-hadrons (easier than by PDG ID)"};

◆ m_motherIDDecoratorKey

SG::WriteDecorHandleKey<xAOD::TruthParticleContainer> DerivationFramework::TruthDecayCollectionMaker::m_motherIDDecoratorKey {this, "motherID", m_outputParticlesKey, "motherID","Name of the decoration which records the ID of the particle's mother"}
private

Definition at line 72 of file TruthDecayCollectionMaker.h.

73{this, "motherID", m_outputParticlesKey, "motherID","Name of the decoration which records the ID of the particle's mother"};

◆ m_originAccessorKey

SG::ReadDecorHandleKey<xAOD::TruthParticleContainer> DerivationFramework::TruthDecayCollectionMaker::m_originAccessorKey {this, "Input_classifierParticleOrigin", m_particlesKey, "classifierParticleOrigin","Name of the decoration which records the particle origin as determined by the MCTruthClassifier"}
private

Definition at line 49 of file TruthDecayCollectionMaker.h.

50{this, "Input_classifierParticleOrigin", m_particlesKey, "classifierParticleOrigin","Name of the decoration which records the particle origin as determined by the MCTruthClassifier"};

◆ m_originDecoratorKey

SG::WriteDecorHandleKey<xAOD::TruthParticleContainer> DerivationFramework::TruthDecayCollectionMaker::m_originDecoratorKey {this, "classifierParticleOrigin", m_outputParticlesKey, "classifierParticleOrigin","Name of the decoration which records the particle origin as determined by the MCTruthClassifier"}
private

FIXME Using WriteDecorHandles for decorations on a Container created in the current algorithm is unnecessary.

Definition at line 64 of file TruthDecayCollectionMaker.h.

65{this, "classifierParticleOrigin", m_outputParticlesKey, "classifierParticleOrigin","Name of the decoration which records the particle origin as determined by the MCTruthClassifier"};

◆ m_outcomeAccessorKey

SG::ReadDecorHandleKey<xAOD::TruthParticleContainer> DerivationFramework::TruthDecayCollectionMaker::m_outcomeAccessorKey {this, "Input_classifierParticleOutCome", m_particlesKey, "classifierParticleOutCome","Name of the decoration which records the particle outcome as determined by the MCTruthClassifier"}
private

Definition at line 53 of file TruthDecayCollectionMaker.h.

54{this, "Input_classifierParticleOutCome", m_particlesKey, "classifierParticleOutCome","Name of the decoration which records the particle outcome as determined by the MCTruthClassifier"};

◆ m_outcomeDecoratorKey

SG::WriteDecorHandleKey<xAOD::TruthParticleContainer> DerivationFramework::TruthDecayCollectionMaker::m_outcomeDecoratorKey {this, "classifierParticleOutCome", m_outputParticlesKey, "classifierParticleOutCome","Name of the decoration which records the particle outcome as determined by the MCTruthClassifier"}
private

Definition at line 68 of file TruthDecayCollectionMaker.h.

69{this, "classifierParticleOutCome", m_outputParticlesKey, "classifierParticleOutCome","Name of the decoration which records the particle outcome as determined by the MCTruthClassifier"};

◆ m_outputParticlesKey

SG::WriteHandleKey<xAOD::TruthParticleContainer> DerivationFramework::TruthDecayCollectionMaker::m_outputParticlesKey {this, "NewParticleKey", "", "WriteHandleKey for new TruthParticleContainer"}
private

Definition at line 61 of file TruthDecayCollectionMaker.h.

62{this, "NewParticleKey", "", "WriteHandleKey for new TruthParticleContainer"};

◆ m_outputVerticesKey

SG::WriteHandleKey<xAOD::TruthVertexContainer> DerivationFramework::TruthDecayCollectionMaker::m_outputVerticesKey {this, "NewVertexKey", "", "WriteHandleKey for new TruthVertexContainer"}
private

Definition at line 59 of file TruthDecayCollectionMaker.h.

60{this, "NewVertexKey", "", "WriteHandleKey for new TruthVertexContainer"};

◆ m_particlesKey

SG::ReadHandleKey<xAOD::TruthParticleContainer> DerivationFramework::TruthDecayCollectionMaker::m_particlesKey {this, "ParticlesKey", "TruthParticles", "ReadHandleKey for input TruthParticleContainer"}
private

Definition at line 47 of file TruthDecayCollectionMaker.h.

48{this, "ParticlesKey", "TruthParticles", "ReadHandleKey for input TruthParticleContainer"};

◆ m_pdgIdsToKeep

Gaudi::Property<std::vector<int> > DerivationFramework::TruthDecayCollectionMaker::m_pdgIdsToKeep {this, "PDGIDsToKeep", {}, "PDG IDs of particles to build the collection from"}
private

< List of PDG IDs to build this collection from

Definition at line 35 of file TruthDecayCollectionMaker.h.

36{this, "PDGIDsToKeep", {}, "PDG IDs of particles to build the collection from"};

◆ m_rejectHadronChildren

Gaudi::Property<bool> DerivationFramework::TruthDecayCollectionMaker::m_rejectHadronChildren {this, "RejectHadronChildren", false, "Drop hadron descendants"}
private

< Option to reject hadron descendants

Definition at line 43 of file TruthDecayCollectionMaker.h.

44{this, "RejectHadronChildren", false, "Drop hadron descendants"};

◆ m_typeAccessorKey

SG::ReadDecorHandleKey<xAOD::TruthParticleContainer> DerivationFramework::TruthDecayCollectionMaker::m_typeAccessorKey {this, "Input_classifierParticleType", m_particlesKey, "classifierParticleType","Name of the decoration which records the particle type as determined by the MCTruthClassifier"}
private

Definition at line 51 of file TruthDecayCollectionMaker.h.

52{this, "Input_classifierParticleType", m_particlesKey, "classifierParticleType","Name of the decoration which records the particle type as determined by the MCTruthClassifier"};

◆ m_typeDecoratorKey

SG::WriteDecorHandleKey<xAOD::TruthParticleContainer> DerivationFramework::TruthDecayCollectionMaker::m_typeDecoratorKey {this, "classifierParticleType", m_outputParticlesKey, "classifierParticleType","Name of the decoration which records the particle type as determined by the MCTruthClassifier"}
private

Definition at line 66 of file TruthDecayCollectionMaker.h.

67{this, "classifierParticleType", m_outputParticlesKey, "classifierParticleType","Name of the decoration which records the particle type as determined by the MCTruthClassifier"};

The documentation for this class was generated from the following files: