ATLAS Offline Software
TruthVertex_v1.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3 */
4 
5 // System include(s):
6 #include <cmath>
7 
8 // xAOD include(s):
10 
11 // Local include(s):
14 
15 namespace xAOD {
16 
18  : SG::AuxElement() {
19 
20  }
21 
23  //
24  // Implementation for the "MC specific" functions
25  //
26 
28  setStatus )
30  setUid )
31 
32  //
34 
36  //
37  // Implementation for the links to the truth particles
38  //
39 
41  incomingParticleLinks,
42  setIncomingParticleLinks )
43 
45  static const SG::AuxElement::Accessor< TruthVertex_v1::TPLinks_t >
46  incomingParticleLinksAcc( "incomingParticleLinks" );
47 
48  size_t TruthVertex_v1::nIncomingParticles() const {
49 
50  // Check if the variable is available:
51  if( ! incomingParticleLinksAcc.isAvailable( *this ) ) {
52  // If not, just tell the user that there aren't any incoming particles:
53  return 0;
54  }
55 
56  // Return the size of the vector:
57  return incomingParticleLinksAcc( *this ).size();
58  }
59 
60  std::vector<const TruthParticle*> TruthVertex_v1::particles_in() const {
61  std::vector<const TruthParticle*> res;
62  for (size_t i=0; i<nIncomingParticles();i++) res.push_back(incomingParticle(i));
63  return res;
64  }
65  std::vector<const TruthParticle*> TruthVertex_v1::particles_out() const {
66  std::vector<const TruthParticle*> res;
67  for (size_t i=0; i<nOutgoingParticles();i++) res.push_back(outgoingParticle(i));
68  return res;
69  }
71 
72  // Check that the variable exists, and that it has enough elements in it:
73  if( ( ! incomingParticleLinksAcc.isAvailable( *this ) ) ||
74  ( incomingParticleLinksAcc( *this ).size() <= index ) ) {
75  return nullptr;
76  }
77 
78  // Retrieve the link object and check its validity:
79  const TPLink_t& ipl = incomingParticleLinksAcc( *this )[ index ];
80  if( ! ipl.isValid() ) {
81  return nullptr;
82  }
83 
84  // Finally, de-reference the link:
85  return *ipl;
86  }
87 
89 
90  incomingParticleLinksAcc( *this ).push_back( link );
91  return;
92  }
93 
95 
96  incomingParticleLinksAcc( *this ).clear();
97  return;
98  }
99 
101  outgoingParticleLinks,
102  setOutgoingParticleLinks )
103 
104 
106  outgoingParticleLinksAcc( "outgoingParticleLinks" );
107 
108  size_t TruthVertex_v1::nOutgoingParticles() const {
109 
110  // Check if the variable is available:
111  if( ! outgoingParticleLinksAcc.isAvailable( *this ) ) {
112  // If not, just tell the user that there aren't any outgoing particles:
113  return 0;
114  }
115 
116  // Return the size of the vector:
117  return outgoingParticleLinksAcc( *this ).size();
118  }
119 
121 
122  // Check that the variable exists, and that it has enough elements in it:
123  if( ( ! outgoingParticleLinksAcc.isAvailable( *this ) ) ||
124  ( outgoingParticleLinksAcc( *this ).size() <= index ) ) {
125  return nullptr;
126  }
127 
128  // Retrieve the link object and check its validity:
129  const TPLink_t& opl = outgoingParticleLinksAcc( *this )[ index ];
130  if( ! opl.isValid() ) {
131  return nullptr;
132  }
133 
134  // Finally, de-reference the link:
135  return *opl;
136  }
137 
139 
140  outgoingParticleLinksAcc( *this ).push_back( link );
141  return;
142  }
143 
145 
146  outgoingParticleLinksAcc( *this ).clear();
147  return;
148  }
149 
150  //
152 
154  //
155  // Implementation of the functions specifying the vertex's position
156  //
157 
159 
161 
163 
165 
166  // Do the calculation by hand. Could make it faster than this even in a
167  // future iteration...
168  return std::sqrt( x() * x() + y() * y() );
169  }
170 
171  float TruthVertex_v1::eta() const {
172 
173  // This is not necessarily what Andy was thinking about...
174  return genvecV4().Eta();
175  }
176 
177  float TruthVertex_v1::phi() const {
178 
179  // This is not necessarily what Andy was thinking about...
180  return genvecV4().Phi();
181  }
182 
184 
186  return FourVec_t(x(), y(), z(), t());
187  }
188 
190  return GenVecFourVec_t(x(), y(), z(), t());
191  }
192 
193  //
195 
197 
198  return Type::TruthVertex;
199  }
200 
202 
203  // Prepare the incoming particle links for persistification:
204  if( incomingParticleLinksAcc.isAvailableWritable( *this ) ) {
205  TPLinks_t::iterator itr = incomingParticleLinksAcc( *this ).begin();
206  TPLinks_t::iterator end = incomingParticleLinksAcc( *this ).end();
207  for( ; itr != end; ++itr ) {
208  itr->toPersistent();
209  }
210  }
211 
212  // Prepare the outgoing particle links for persistification:
213  if( outgoingParticleLinksAcc.isAvailableWritable( *this ) ) {
214  TPLinks_t::iterator itr = outgoingParticleLinksAcc( *this ).begin();
215  TPLinks_t::iterator end = outgoingParticleLinksAcc( *this ).end();
216  for( ; itr != end; ++itr ) {
217  itr->toPersistent();
218  }
219  }
220 
221  return;
222  }
223 
224 } // namespace xAOD
xAOD::iterator
JetConstituentVector::iterator iterator
Definition: JetConstituentVector.cxx:68
xAOD::TruthVertex_v1::clearOutgoingParticleLinks
void clearOutgoingParticleLinks()
Remove all outgoing particles.
Definition: TruthVertex_v1.cxx:144
xAOD::TruthVertex_v1::nOutgoingParticles
size_t nOutgoingParticles() const
Get the number of outgoing particles.
xAOD::AUXSTORE_PRIMITIVE_SETTER_AND_GETTER
AUXSTORE_PRIMITIVE_SETTER_AND_GETTER(BTagging_v1, float, IP2D_pb, setIP2D_pb) AUXSTORE_PRIMITIVE_SETTER_AND_GETTER(BTagging_v1
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
SG
Forward declaration.
Definition: CaloCellPacker_400_500.h:32
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
SG::Accessor
Helper class to provide type-safe access to aux data.
Definition: Control/AthContainers/AthContainers/Accessor.h:68
perp
Scalar perp() const
perp method - perpenticular length
Definition: AmgMatrixBasePlugin.h:44
AuxStoreAccessorMacros.h
xAOD::TruthVertex_v1::clearIncomingParticleLinks
void clearIncomingParticleLinks()
Remove all incoming particles.
Definition: TruthVertex_v1.cxx:94
index
Definition: index.py:1
TruthParticleContainer.h
xAOD
ICaloAffectedTool is abstract interface for tools checking if 4 mom is in calo affected region.
Definition: ICaloAffectedTool.h:24
TruthVertex_v1.h
read_hist_ntuple.t
t
Definition: read_hist_ntuple.py:5
xAOD::TruthVertex_v1::addOutgoingParticleLink
void addOutgoingParticleLink(const TPLink_t &link)
Add one outgoing particle.
Definition: TruthVertex_v1.cxx:138
const
bool const RAWDATA *ch2 const
Definition: LArRodBlockPhysicsV0.cxx:560
x
#define x
xAOD::TruthVertex_v1::y
float y() const
Vertex y displacement.
xAOD::TruthVertex_v1::t
float t() const
Vertex time.
mergePhysValFiles.end
end
Definition: DataQuality/DataQualityUtils/scripts/mergePhysValFiles.py:92
JetVar::Accessor
SG::AuxElement::Accessor< T > Accessor
Definition: JetVariable.h:31
SG::IAuxElement::index
size_t index() const
Return the index of this element within its container.
xAOD::uid
uid
Definition: TruthVertex_v1.cxx:29
lumiFormat.i
int i
Definition: lumiFormat.py:85
z
#define z
xAOD::float
float
Definition: BTagging_v1.cxx:168
xAOD::TruthParticle_v1
Class describing a truth particle in the MC record.
Definition: TruthParticle_v1.h:37
res
std::pair< std::vector< unsigned int >, bool > res
Definition: JetGroupProductTest.cxx:11
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
TruthVertex_v1 TruthVertex
Typedef to implementation.
Definition: TruthVertex.h:15
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
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.
y
#define y
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.
merge.status
status
Definition: merge.py:16
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::AUXSTORE_OBJECT_SETTER_AND_GETTER
AUXSTORE_OBJECT_SETTER_AND_GETTER(CaloRings_v1, RingSetLinks, ringSetLinks, setRingSetLinks) unsigned CaloRings_v1
Definition: CaloRings_v1.cxx:27