7#ifndef XAODTRUTH_VERSIONS_TRUTHVERTEX_V1_H
8#define XAODTRUTH_VERSIONS_TRUTHVERTEX_V1_H
11#include <TLorentzVector.h>
15#include "AthLinks/ElementLink.h"
24#include "Math/Vector4D.h"
143 Type::ObjectType
type()
const;
150 if (!v) { os <<
"Vtx: Empty vertex" << std::endl;
return os;}
152 os << v->uid() <<
" status=";
154 os <<
" (x,y,z,t)=" << v->x() <<
"," << v->y() <<
"," << v->z() <<
"," << v->t();
Base class for elements of a container that can have aux data.
#define SG_BASE(D, B)
Declare that class D derives from class B.
Base class for elements of a container that can have aux data.
Class describing a truth particle in the MC record.
Class describing a truth vertex in the MC record.
const TruthParticle_v1 * outgoingParticle(size_t index) const
Get one of the outgoing particles.
const TPLinks_t & outgoingParticleLinks() const
Get all 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.
float eta() const
Vertex pseudorapidity.
FourVec_t v4() const
The full 4-vector of the vertex.
void setStatus(int value)
Set the vertex status.
int status() const
Get the vertex status.
void setZ(float value)
Set the vertex's longitudinal distance from the origin.
float y() const
Vertex y displacement.
void clearIncomingParticleLinks()
Remove all incoming particles.
TruthVertex_v1()
Default constructor.
float phi() const
Vertex azimuthal angle.
Type::ObjectType type() const
The type of the object as a simple enumeration.
int uid() const
Get the vertex unique ID.
void setUid(int value)
Set the vertex unique ID.
void setOutgoingParticleLinks(const TPLinks_t &links)
Set all the outgoing particles.
ElementLink< TruthParticleContainer > TPLink_t
Type of one truth particle link.
const TruthParticle_v1 * incomingParticle(size_t index) const
Get one of the incoming particles.
void setIncomingParticleLinks(const TPLinks_t &links)
Set all the incoming particles.
float t() const
Vertex time.
std::vector< const TruthParticle * > particles_out() const
Get the outgoing particles.
float perp() const
Vertex transverse distance from the beam line.
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > GenVecFourVec_t
Base 4 Momentum type for TruthVector.
size_t nOutgoingParticles() const
Get the number of outgoing particles.
TLorentzVector FourVec_t
The 4-vector type.
size_t nIncomingParticles() const
Get the number of incoming particles.
void toPersistent()
Function making sure that the object is ready for persistification.
void setT(float value)
Set the vertex time.
float x() const
Vertex x displacement.
std::vector< const TruthParticle * > particles_in() const
Get the incoming particles.
void setX(float value)
Set the x displacement of the vertex.
void addIncomingParticleLink(const TPLink_t &link)
Add one incoming particle.
void clearOutgoingParticleLinks()
Remove all outgoing particles.
std::vector< TPLink_t > TPLinks_t
Type used to save the links to incoming and outgoing particles.
void setY(float value)
Set the y displacement of the vertex.
const TPLinks_t & incomingParticleLinks() const
Get all the incoming particles.
GenVecFourVec_t genvecV4() const
The full 4-vector of the particle : GenVector form.
ICaloAffectedTool is abstract interface for tools checking if 4 mom is in calo affected region.
std::ostream & operator<<(std::ostream &out, const std::pair< FIRST, SECOND > &pair)
Helper print operator.