ATLAS Offline Software
Loading...
Searching...
No Matches
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
26namespace xAOD {
27
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
69 void setIncomingParticleLinks( const TPLinks_t& links );
71 size_t nIncomingParticles() const;
73 const TruthParticle_v1* incomingParticle( size_t index ) const;
75 void addIncomingParticleLink( const TPLink_t& link );
78
82 void setOutgoingParticleLinks( const TPLinks_t& links );
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
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
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.
Definition AuxElement.h:483
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.
Definition index.py:1
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.