ATLAS Offline Software
TruthVertex_v1.h
Go to the documentation of this file.
1 // -*- C++ -*-
2 
3 /*
4  Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
5 */
6 
7 #ifndef XAODTRUTH_VERSIONS_TRUTHVERTEX_V1_H
8 #define XAODTRUTH_VERSIONS_TRUTHVERTEX_V1_H
9 
10 // ROOT include(s):
11 #include <TLorentzVector.h>
12 
13 // EDM include(s):
15 #include "AthLinks/ElementLink.h"
16 
17 // xAOD include(s):
18 #include "xAODBase/ObjectType.h"
19 
20 // Local include(s):
22 
23 // ROOT include(s):
24 #include "Math/Vector4D.h"
25 
26 namespace xAOD {
27 
37  class TruthVertex_v1 : public SG::AuxElement {
38 
39  public:
42 
45 
47  int status() const;
49  void setStatus( int value );
50 
52  int uid() const;
54  void setUid( int value );
55 
57 
60 
64  typedef std::vector< TPLink_t > TPLinks_t;
65 
71  size_t nIncomingParticles() const;
73  const TruthParticle_v1* incomingParticle( size_t index ) const;
75  void addIncomingParticleLink( const TPLink_t& link );
78 
84  size_t nOutgoingParticles() const;
86  const TruthParticle_v1* outgoingParticle( size_t index ) const;
88  std::vector<const TruthParticle*> particles_in() const;
90  std::vector<const TruthParticle*> particles_out() const;
92  void addOutgoingParticleLink( const TPLink_t& link );
95 
97 
100 
102  float x() const;
104  void setX( float value );
105 
107  float y() const;
109  void setY( float value );
110 
112  float z() const;
114  void setZ( float value );
115 
117  float perp() const;
119  float eta() const;
121  float phi() const;
122 
124  float t() const;
126  void setT( float value );
127 
129  typedef TLorentzVector FourVec_t;
130 
132  FourVec_t v4() const;
133 
135  typedef ROOT::Math::LorentzVector<ROOT::Math::PxPyPzE4D<double> > GenVecFourVec_t;
136 
138  GenVecFourVec_t genvecV4() const;
139 
141 
143  Type::ObjectType type() const;
144 
146  void toPersistent();
147 
148  }; // class TruthVertex_v1
149  inline std::ostream& operator<<(std::ostream& os, const TruthVertex_v1* v) {
150  if (!v) { os << "Vtx: Empty vertex" << std::endl; return os;}
151  os << "Vtx: uid=";
152  os << v->uid() << " status=";
153  os << v->status();
154  os << " (x,y,z,t)=" << v->x() << "," << v->y() << "," << v->z() << "," << v->t();
155  /* os << std::endl; AV:Not clear if we need a new line here */
156  return os;
157  }
158 } // namespace xAOD
159 
160 // Declare the inheritance of the type to StoreGate:
161 #include "xAODCore/BaseInfo.h"
163 
164 #endif // XAODTRUTH_VERSIONS_TRUTHVERTEX_V1_H
xAOD::TruthVertex_v1::clearOutgoingParticleLinks
void clearOutgoingParticleLinks()
Remove all outgoing particles.
Definition: TruthVertex_v1.cxx:144
TruthParticleContainerFwd.h
xAOD::TruthVertex_v1::nOutgoingParticles
size_t nOutgoingParticles() const
Get the number of outgoing particles.
xAOD::TruthVertex_v1::GenVecFourVec_t
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > GenVecFourVec_t
Base 4 Momentum type for TruthVector.
Definition: TruthVertex_v1.h:135
xAOD::TruthVertex_v1::phi
float phi() const
Vertex azimuthal angle.
Definition: TruthVertex_v1.cxx:177
xAOD::TruthVertex_v1::FourVec_t
TLorentzVector FourVec_t
The 4-vector type.
Definition: TruthVertex_v1.h:129
xAOD::TruthVertex_v1::particles_in
std::vector< const TruthParticle * > particles_in() const
Get the incoming particles.
Definition: TruthVertex_v1.cxx:60
xAOD::TruthVertex_v1::clearIncomingParticleLinks
void clearIncomingParticleLinks()
Remove all incoming particles.
Definition: TruthVertex_v1.cxx:94
index
Definition: index.py:1
SG::AuxElement
Base class for elements of a container that can have aux data.
Definition: AuxElement.h:483
BaseInfo.h
xAOD::TruthVertex_v1::setT
void setT(float value)
Set the vertex time.
athena.value
value
Definition: athena.py:124
xAOD
ICaloAffectedTool is abstract interface for tools checking if 4 mom is in calo affected region.
Definition: ICaloAffectedTool.h:24
xAOD::TruthVertex_v1::setIncomingParticleLinks
void setIncomingParticleLinks(const TPLinks_t &links)
Set all the incoming particles.
xAOD::TruthVertex_v1::addOutgoingParticleLink
void addOutgoingParticleLink(const TPLink_t &link)
Add one outgoing particle.
Definition: TruthVertex_v1.cxx:138
xAOD::TruthVertex_v1::setY
void setY(float value)
Set the y displacement of the vertex.
xAOD::TruthVertex_v1::uid
int uid() const
Get the vertex unique ID.
SG_BASE
SG_BASE(xAOD::TruthVertex_v1, SG::AuxElement)
xAOD::TruthVertex_v1::setOutgoingParticleLinks
void setOutgoingParticleLinks(const TPLinks_t &links)
Set all the outgoing particles.
xAOD::TruthVertex_v1::outgoingParticleLinks
const TPLinks_t & outgoingParticleLinks() const
Get all the outgoing particles.
xAOD::TruthVertex_v1::y
float y() const
Vertex y displacement.
xAOD::TruthVertex_v1::t
float t() const
Vertex time.
xAOD::TruthVertex_v1::v4
FourVec_t v4() const
The full 4-vector of the vertex.
Definition: TruthVertex_v1.cxx:185
ObjectType.h
DMTest::links
links
Definition: CLinks_v1.cxx:22
xAOD::TruthParticle_v1
Class describing a truth particle in the MC record.
Definition: TruthParticle_v1.h:37
xAOD::TruthVertex_v1::perp
float perp() const
Vertex transverse distance from the beam line.
Definition: TruthVertex_v1.cxx:164
xAOD::TruthVertex_v1::incomingParticle
const TruthParticle_v1 * incomingParticle(size_t index) const
Get one of the incoming particles.
Definition: TruthVertex_v1.cxx:70
xAOD::TruthVertex_v1::toPersistent
void toPersistent()
Function making sure that the object is ready for persistification.
Definition: TruthVertex_v1.cxx:201
xAOD::TruthVertex_v1::TPLinks_t
std::vector< TPLink_t > TPLinks_t
Type used to save the links to incoming and outgoing particles.
Definition: TruthVertex_v1.h:64
ReadFromCoolCompare.os
os
Definition: ReadFromCoolCompare.py:231
xAOD::TruthVertex_v1
Class describing a truth vertex in the MC record.
Definition: TruthVertex_v1.h:37
xAOD::TruthVertex_v1::particles_out
std::vector< const TruthParticle * > particles_out() const
Get the outgoing particles.
Definition: TruthVertex_v1.cxx:65
xAOD::TruthVertex_v1::x
float x() const
Vertex x displacement.
python.PyAthena.v
v
Definition: PyAthena.py:154
xAOD::TruthVertex_v1::setZ
void setZ(float value)
Set the vertex's longitudinal distance from the origin.
xAOD::TruthVertex_v1::status
int status() const
Get the vertex status.
xAOD::TruthVertex_v1::TPLink_t
ElementLink< TruthParticleContainer > TPLink_t
Type of one truth particle link.
Definition: TruthVertex_v1.h:62
xAOD::TruthVertex_v1::setUid
void setUid(int value)
Set the vertex unique ID.
xAOD::TruthVertex_v1::eta
float eta() const
Vertex pseudorapidity.
Definition: TruthVertex_v1.cxx:171
xAOD::TruthVertex_v1::nIncomingParticles
size_t nIncomingParticles() const
Get the number of incoming particles.
xAOD::TruthVertex_v1::type
Type::ObjectType type() const
The type of the object as a simple enumeration.
Definition: TruthVertex_v1.cxx:196
xAOD::TruthVertex_v1::TruthVertex_v1
TruthVertex_v1()
Default constructor.
Definition: TruthVertex_v1.cxx:17
xAOD::TruthVertex_v1::z
float z() const
Vertex longitudinal distance along the beam line form the origin.
xAOD::TruthVertex_v1::genvecV4
GenVecFourVec_t genvecV4() const
The full 4-vector of the particle : GenVector form.
Definition: TruthVertex_v1.cxx:189
xAODType::ObjectType
ObjectType
Type of objects that have a representation in the xAOD EDM.
Definition: ObjectType.h:32
xAOD::TruthVertex_v1::addIncomingParticleLink
void addIncomingParticleLink(const TPLink_t &link)
Add one incoming particle.
Definition: TruthVertex_v1.cxx:88
xAOD::TruthVertex_v1::outgoingParticle
const TruthParticle_v1 * outgoingParticle(size_t index) const
Get one of the outgoing particles.
Definition: TruthVertex_v1.cxx:120
xAOD::TruthVertex_v1::setStatus
void setStatus(int value)
Set the vertex status.
AuxElement.h
Base class for elements of a container that can have aux data.
xAOD::operator<<
std::ostream & operator<<(std::ostream &out, const std::pair< FIRST, SECOND > &pair)
Helper print operator.
Definition: RDataSource.cxx:53
xAOD::TruthVertex_v1::setX
void setX(float value)
Set the x displacement of the vertex.
xAOD::TruthVertex_v1::incomingParticleLinks
const TPLinks_t & incomingParticleLinks() const
Get all the incoming particles.