|
ATLAS Offline Software
|
Go to the documentation of this file.
8 #ifndef XAODTRUTH_VERSIONS_TRUTHVERTEX_V1_H
9 #define XAODTRUTH_VERSIONS_TRUTHVERTEX_V1_H
12 #include <TLorentzVector.h>
16 #include "AthLinks/ElementLink.h"
25 #include "Math/Vector4D.h"
170 #endif // XAODTRUTH_VERSIONS_TRUTHVERTEX_V1_H
void clearOutgoingParticleLinks()
Remove all outgoing particles.
size_t nOutgoingParticles() const
Get the number of outgoing particles.
int id() const
Obsolete function Vertex ID code HepMC2 id == HepMC3 status, i.e.
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > GenVecFourVec_t
Base 4 Momentum type for TruthVector.
float phi() const
Vertex azimuthal angle.
TLorentzVector FourVec_t
The 4-vector type.
std::vector< const TruthParticle * > particles_in() const
Get the incoming particles.
void clearIncomingParticleLinks()
Remove all incoming particles.
Base class for elements of a container that can have aux data.
void setT(float value)
Set the vertex time.
ICaloAffectedTool is abstract interface for tools checking if 4 mom is in calo affected region.
void setIncomingParticleLinks(const TPLinks_t &links)
Set all the incoming particles.
void addOutgoingParticleLink(const TPLink_t &link)
Add one outgoing particle.
void setY(float value)
Set the y displacement of the vertex.
SG_BASE(xAOD::TruthVertex_v1, SG::AuxElement)
void setOutgoingParticleLinks(const TPLinks_t &links)
Set all the outgoing particles.
const TPLinks_t & outgoingParticleLinks() const
Get all the outgoing particles.
float y() const
Vertex y displacement.
float t() const
Vertex time.
FourVec_t v4() const
The full 4-vector of the vertex.
int barcode() const
Barcode.
Class describing a truth particle in the MC record.
float perp() const
Vertex transverse distance from the beam line.
const TruthParticle_v1 * incomingParticle(size_t index) const
Get one of the incoming particles.
void toPersistent()
Function making sure that the object is ready for persistification.
std::vector< TPLink_t > TPLinks_t
Type used to save the links to incoming and outgoing particles.
Class describing a truth vertex in the MC record.
std::vector< const TruthParticle * > particles_out() const
Get the outgoing particles.
void setId(int value)
Obsolete function Set vertex ID code HepMC2 id == HepMC3 status, i.e.
float x() const
Vertex x displacement.
void setZ(float value)
Set the vertex's longitudinal distance from the origin.
int status() const
New function Vertex status HepMC2 id == HepMC3 status, i.e.
ElementLink< TruthParticleContainer > TPLink_t
Type of one truth particle link.
void setBarcode(int value)
Set barcode.
float eta() const
Vertex pseudorapidity.
size_t nIncomingParticles() const
Get the number of incoming particles.
Type::ObjectType type() const
The type of the object as a simple enumeration.
TruthVertex_v1()
Default constructor.
float z() const
Vertex longitudinal distance along the beam line form the origin.
GenVecFourVec_t genvecV4() const
The full 4-vector of the particle : GenVector form.
ObjectType
Type of objects that have a representation in the xAOD EDM.
void addIncomingParticleLink(const TPLink_t &link)
Add one incoming particle.
const TruthParticle_v1 * outgoingParticle(size_t index) const
Get one of the outgoing particles.
void setStatus(int value)
New function Set vertex status HepMC2 id == HepMC3 status, i.e.
Base class for elements of a container that can have aux data.
void setX(float value)
Set the x displacement of the vertex.
const TPLinks_t & incomingParticleLinks() const
Get all the incoming particles.