|
ATLAS Offline Software
|
#include <CollectionMakerHelpers.h>
|
static int | addTruthParticle (const xAOD::TruthParticle &oldPart, xAOD::TruthParticleContainer *partCont, xAOD::TruthVertexContainer *vertCont, std::vector< int > &seenParticles, const int generations, bool includeVertex=true) |
| < Helper function to add truth particles to a collection More...
|
|
static int | addTruthVertex (const xAOD::TruthVertex &oldVert, xAOD::TruthParticleContainer *partCont, xAOD::TruthVertexContainer *vertCont, std::vector< int > &seenParticles, const int generations) |
| Helper function to set up a truth particle based on an old particle. More...
|
|
static xAOD::TruthParticle * | setupTruthParticle (const xAOD::TruthParticle &oldPart, xAOD::TruthParticleContainer *partCont) |
|
Definition at line 20 of file CollectionMakerHelpers.h.
◆ addTruthParticle()
< Helper function to add truth particles to a collection
Helper function to add truth vertices to a collection
Definition at line 55 of file CollectionMakerHelpers.cxx.
60 for (
size_t p=0;
p<partCont->
size();++
p){
71 int myIndex = partCont->
size()-1;
78 (*vertCont)[vertIndex]->addIncomingParticleLink( eltp );
◆ addTruthVertex()
Helper function to set up a truth particle based on an old particle.
Definition at line 25 of file CollectionMakerHelpers.cxx.
32 int myIndex = vertCont->
size()-1;
37 xTruthVertex->
setX(oldVert.
x());
38 xTruthVertex->
setY(oldVert.
y());
39 xTruthVertex->
setZ(oldVert.
z());
40 xTruthVertex->
setT(oldVert.
t());
49 (*partCont)[partIndex]->setProdVtxLink( eltv );
◆ setupTruthParticle()
Definition at line 84 of file CollectionMakerHelpers.cxx.
103 xTruthParticle->
setM(oldPart.
m());
104 xTruthParticle->
setPx(oldPart.
px());
105 xTruthParticle->
setPy(oldPart.
py());
106 xTruthParticle->
setPz(oldPart.
pz());
107 xTruthParticle->
setE(oldPart.
e());
114 typeDecorator(*xTruthParticle) = classifierParticleTypeAcc.withDefault(oldPart, 0);
115 originDecorator(*xTruthParticle) = classifierParticleOriginAcc.withDefault(oldPart, 0);
116 outcomeDecorator(*xTruthParticle) = classifierParticleOutComeAcc.withDefault(oldPart, 0);
117 return xTruthParticle;
The documentation for this class was generated from the following files:
size_t nOutgoingParticles() const
Get the number of outgoing particles.
void setStatus(int value)
Set status code.
static int addTruthParticle(const xAOD::TruthParticle &oldPart, xAOD::TruthParticleContainer *partCont, xAOD::TruthVertexContainer *vertCont, std::vector< int > &seenParticles, const int generations, bool includeVertex=true)
< Helper function to add truth particles to a collection
std::string find(const std::string &s)
return a remapped string
float pz() const
The z component of the particle's momentum.
void setE(float value)
Set the energy of the particle.
void setBarcode(int value)
Set barcode.
bool setPolarizationParameter(float value, PolParam parameter)
Set method for polarization parameter values.
bool valid() const
Check if the stored values are valid.
Polarization polarization() const
Retrieve a full Polarization with a single call.
float px() const
The x component of the particle's momentum.
void setT(float value)
Set the vertex time.
void setPx(float value)
Set the x component of the particle's momentum.
float py() const
The y component of the particle's momentum.
void addOutgoingParticleLink(const TPLink_t &link)
Add one outgoing particle.
void setY(float value)
Set the y displacement of the vertex.
bool is_same_particle(const T1 &p1, const T2 &p2)
Method to establish if two particles in the GenEvent actually represent the same particle.
float y() const
Vertex y displacement.
int generations(const T &p)
Method to return how many interactions a particle has undergone during simulation (only to be used in...
float t() const
Vertex time.
bool hasDecayVtx() const
Check for a decay vertex on this particle.
void setM(float value)
Also store the mass.
virtual double e() const override final
The total energy of the particle.
Helper class to provide type-safe access to aux data.
Class describing a truth particle in the MC record.
TruthParticle_v1 TruthParticle
Typedef to implementation.
void setPy(float value)
Set the y component of the particle's momentum.
TruthVertex_v1 TruthVertex
Typedef to implementation.
ElementLink implementation for ROOT usage.
const TruthVertex_v1 * decayVtx() const
The decay vertex of this particle.
Class describing a truth vertex in the MC record.
void setPdgId(int pid)
Set PDG ID code.
value_type push_back(value_type pElem)
Add an element to the end of the collection.
void setId(int value)
Obsolete function Set vertex ID code HepMC2 id == HepMC3 status, i.e.
bool polarizationParameter(float &value, PolParam parameter) const
Accessor for polarization parameters.
float x() const
Vertex x displacement.
void setZ(float value)
Set the vertex's longitudinal distance from the origin.
int status() const
Status code.
void setBarcode(int value)
Set barcode.
@ polarizationTheta
Polarization in ( )
@ polarizationPhi
Polarization in ( )
float z() const
Vertex longitudinal distance along the beam line form the origin.
void setPz(float value)
Set the z component of the particle's momentum.
void setDecayVtxLink(const ElementLink< TruthVertexContainer > &link)
Set the decay vertex of the particle.
int pdgId() const
PDG ID code.
size_type size() const noexcept
Returns the number of elements in the collection.
const TruthParticle_v1 * outgoingParticle(size_t index) const
Get one of the outgoing particles.
static xAOD::TruthParticle * setupTruthParticle(const xAOD::TruthParticle &oldPart, xAOD::TruthParticleContainer *partCont)
virtual double m() const override final
The mass of the particle.
static int addTruthVertex(const xAOD::TruthVertex &oldVert, xAOD::TruthParticleContainer *partCont, xAOD::TruthVertexContainer *vertCont, std::vector< int > &seenParticles, const int generations)
Helper function to set up a truth particle based on an old particle.
void setX(float value)
Set the x displacement of the vertex.