ATLAS Offline Software
CollectionMakerHelpers.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3 */
4 
6 // CollectionMakerHelpers.cxx
7 // Helper functions used in a couple of collection makers to add
8 // particles and vertices to light collections
9 
10 #include "CollectionMakerHelpers.h"
17 
19 
20 // STL includes
21 #include <vector>
22 // For a find in the vector
23 #include <algorithm>
24 
26  xAOD::TruthVertexContainer* vertCont, std::vector<int>& seenParticles,
27  const int generations) {
28  // Make a new vertex and add it to the container
29  xAOD::TruthVertex* xTruthVertex = new xAOD::TruthVertex();
30  vertCont->push_back( xTruthVertex );
31  // Get a link to this vertex -- will be used to set production vertices on all the next particles
32  int myIndex = vertCont->size()-1;
33  ElementLink<xAOD::TruthVertexContainer> eltv(*vertCont, myIndex);
34  // Set properties
35  xTruthVertex->setId(HepMC::status(oldVert));
36  xTruthVertex->setBarcode(HepMC::barcode(&oldVert)); // FIXME barcode-based
37  xTruthVertex->setX(oldVert.x());
38  xTruthVertex->setY(oldVert.y());
39  xTruthVertex->setZ(oldVert.z());
40  xTruthVertex->setT(oldVert.t());
41  // If we are done, then stop here
42  if (generations==0) return myIndex;
43  // Add all the outgoing particles
44  for (size_t n=0;n<oldVert.nOutgoingParticles();++n){
45  if (!oldVert.outgoingParticle(n)) continue; // Just in case we removed some truth particles, e.g. G4 decays
46  int partIndex = CollectionMakerHelpers::addTruthParticle( *oldVert.outgoingParticle(n), partCont, vertCont, seenParticles, generations-1);
47  ElementLink<xAOD::TruthParticleContainer> eltp( *partCont, partIndex);
48  xTruthVertex->addOutgoingParticleLink( eltp );
49  (*partCont)[partIndex]->setProdVtxLink( eltv );
50  }
51  // Return a link to this vertex
52  return myIndex;
53 }
54 
56  xAOD::TruthVertexContainer* vertCont, std::vector<int>& seenParticles,
57  const int generations, bool includeVertex) {
58  // See if we've seen it - note, could also do this with a unary function on the container itself
59  if (std::find(seenParticles.begin(), seenParticles.end(), HepMC::uniqueID(&oldPart)) != seenParticles.end()){
60  for (size_t p=0;p<partCont->size();++p){
61  // Was it a hit?
62  const xAOD::TruthParticle *theParticle = (*partCont)[p];
63  if (HepMC::uniqueID(theParticle) == HepMC::uniqueID(&oldPart)) return p;
64  } // Look through the old container
65  } // Found it in the old container
66  // Now we have seen it
67  seenParticles.push_back(HepMC::uniqueID(&oldPart));
68  // Make a nw truth particle
69  xAOD::TruthParticle* xTruthParticle = setupTruthParticle(oldPart,partCont);
70  // Make a link to this particle
71  int myIndex = partCont->size()-1;
72  ElementLink<xAOD::TruthParticleContainer> eltp(*partCont, myIndex);
73  // Decay vertex information
74  if (oldPart.hasDecayVtx() && includeVertex) {
75  int vertIndex = CollectionMakerHelpers::addTruthVertex( *oldPart.decayVtx(), partCont, vertCont, seenParticles, generations);
76  ElementLink<xAOD::TruthVertexContainer> eltv( *vertCont, vertIndex );
77  xTruthParticle->setDecayVtxLink( eltv );
78  (*vertCont)[vertIndex]->addIncomingParticleLink( eltp );
79  }
80  // Return a link to this particle
81  return myIndex;
82 }
83 
85  // Set up decorators
86  const static SG::Decorator< unsigned int > originDecorator("classifierParticleOrigin");
87  const static SG::Decorator< unsigned int > typeDecorator("classifierParticleType");
88  const static SG::Decorator< unsigned int > outcomeDecorator("classifierParticleOutCome");
89  const static SG::Decorator< int > motherIDDecorator("motherID");
90  const static SG::Decorator< int > daughterIDDecorator("daughterID");
91 
92  static const SG::ConstAccessor<unsigned int> classifierParticleTypeAcc("classifierParticleType");
93  static const SG::ConstAccessor<unsigned int> classifierParticleOriginAcc("classifierParticleOrigin");
94  static const SG::ConstAccessor<unsigned int> classifierParticleOutComeAcc("classifierParticleOutCome");
95 
96  // Make a truth particle and add it to the container
97  xAOD::TruthParticle* xTruthParticle = new xAOD::TruthParticle();
98  partCont->push_back( xTruthParticle );
99  // Fill with numerical content
100  xTruthParticle->setPdgId(oldPart.pdgId());
101  xTruthParticle->setBarcode(HepMC::barcode(&oldPart)); // FIXME barcode-based
102  xTruthParticle->setStatus(oldPart.status());
103  xTruthParticle->setM(oldPart.m());
104  xTruthParticle->setPx(oldPart.px());
105  xTruthParticle->setPy(oldPart.py());
106  xTruthParticle->setPz(oldPart.pz());
107  xTruthParticle->setE(oldPart.e());
108  // Copy over the polarization information if it's there
109  if (oldPart.polarization().valid()){
112  }
113  // Copy over the decorations if they are available
114  typeDecorator(*xTruthParticle) = classifierParticleTypeAcc.withDefault(oldPart, 0);
115  originDecorator(*xTruthParticle) = classifierParticleOriginAcc.withDefault(oldPart, 0);
116  outcomeDecorator(*xTruthParticle) = classifierParticleOutComeAcc.withDefault(oldPart, 0);
117  return xTruthParticle;
118 }
CollectionMakerHelpers.h
xAOD::TruthVertex_v1::nOutgoingParticles
size_t nOutgoingParticles() const
Get the number of outgoing particles.
xAOD::TruthParticle_v1::setStatus
void setStatus(int value)
Set status code.
python.PerfMonSerializer.p
def p
Definition: PerfMonSerializer.py:743
DerivationFramework::CollectionMakerHelpers::addTruthParticle
static int addTruthParticle(const xAOD::TruthParticle &oldPart, xAOD::TruthParticleContainer *partCont, xAOD::TruthVertexContainer *vertCont, std::vector< int > &seenParticles, const int generations, bool includeVertex=true)
< Helper function to add truth particles to a collection
Definition: CollectionMakerHelpers.cxx:55
find
std::string find(const std::string &s)
return a remapped string
Definition: hcg.cxx:135
xAOD::TruthParticle_v1::pz
float pz() const
The z component of the particle's momentum.
xAOD::TruthParticle_v1::setE
void setE(float value)
Set the energy of the particle.
Definition: TruthParticle_v1.cxx:235
xAOD::TruthParticle_v1::setBarcode
void setBarcode(int value)
Set barcode.
TruthVertexContainer.h
TruthParticleContainer.h
xAOD::TruthParticle_v1::setPolarizationParameter
bool setPolarizationParameter(float value, PolParam parameter)
Set method for polarization parameter values.
Definition: TruthParticle_v1.cxx:351
xAOD::TruthParticle_v1::Polarization::valid
bool valid() const
Check if the stored values are valid.
Definition: TruthParticle_v1.h:379
xAOD::TruthParticle_v1::polarization
Polarization polarization() const
Retrieve a full Polarization with a single call.
Definition: TruthParticle_v1.cxx:383
xAOD::TruthParticle_v1::px
float px() const
The x component of the particle's momentum.
xAOD::TruthVertex_v1::setT
void setT(float value)
Set the vertex time.
xAOD::TruthParticle_v1::setPx
void setPx(float value)
Set the x component of the particle's momentum.
xAOD::TruthParticle_v1::py
float py() const
The y component of the particle's momentum.
SG::ConstAccessor< unsigned int >
xAOD::TruthVertex_v1::addOutgoingParticleLink
void addOutgoingParticleLink(const TPLink_t &link)
Add one outgoing particle.
Definition: TruthVertex_v1.cxx:139
xAOD::TruthVertex_v1::setY
void setY(float value)
Set the y displacement of the vertex.
xAOD::TruthVertex_v1::y
float y() const
Vertex y displacement.
HepMC::BarcodeBased::generations
int generations(const T &p)
Method to return how many interactions a particle has undergone during simulation (only to be used in...
Definition: MagicNumbers.h:204
xAOD::TruthVertex_v1::t
float t() const
Vertex time.
TruthParticleAuxContainer.h
xAOD::TruthParticle_v1::hasDecayVtx
bool hasDecayVtx() const
Check for a decay vertex on this particle.
xAOD::TruthParticle_v1::setM
void setM(float value)
Also store the mass.
Definition: TruthParticle_v1.cxx:241
xAOD::TruthParticle_v1::e
virtual double e() const override final
The total energy of the particle.
SG::Decorator
Helper class to provide type-safe access to aux data.
Definition: Decorator.h:58
beamspotman.n
n
Definition: beamspotman.py:731
HepMC::barcode
int barcode(const T *p)
Definition: Barcode.h:16
xAOD::TruthParticle_v1
Class describing a truth particle in the MC record.
Definition: TruthParticle_v1.h:41
xAOD::TruthParticle
TruthParticle_v1 TruthParticle
Typedef to implementation.
Definition: Event/xAOD/xAODTruth/xAODTruth/TruthParticle.h:15
HepMC::uniqueID
int uniqueID(const T &p)
Definition: MagicNumbers.h:113
xAOD::TruthParticle_v1::setPy
void setPy(float value)
Set the y component of the particle's momentum.
xAOD::TruthVertex
TruthVertex_v1 TruthVertex
Typedef to implementation.
Definition: TruthVertex.h:15
DataVector
Derived DataVector<T>.
Definition: DataVector.h:581
TruthVertexAuxContainer.h
xAOD::TruthParticle_v1::decayVtx
const TruthVertex_v1 * decayVtx() const
The decay vertex of this particle.
xAOD::TruthVertex_v1
Class describing a truth vertex in the MC record.
Definition: TruthVertex_v1.h:41
xAOD::TruthParticle_v1::setPdgId
void setPdgId(int pid)
Set PDG ID code.
DataVector::push_back
value_type push_back(value_type pElem)
Add an element to the end of the collection.
xAOD::TruthVertex_v1::setId
void setId(int value)
Obsolete function Set vertex ID code HepMC2 id == HepMC3 status, i.e.
xAOD::TruthParticle_v1::polarizationParameter
bool polarizationParameter(float &value, PolParam parameter) const
Accessor for polarization parameters.
Definition: TruthParticle_v1.cxx:328
xAOD::TruthVertex_v1::x
float x() const
Vertex x displacement.
xAOD::TruthVertex_v1::setZ
void setZ(float value)
Set the vertex's longitudinal distance from the origin.
xAOD::TruthParticle_v1::status
int status() const
Status code.
xAOD::TruthVertex_v1::setBarcode
void setBarcode(int value)
Set barcode.
xAOD::TruthParticle_v1::polarizationTheta
@ polarizationTheta
Polarization in ( )
Definition: TruthParticle_v1.h:323
xAOD::TruthParticle_v1::polarizationPhi
@ polarizationPhi
Polarization in ( )
Definition: TruthParticle_v1.h:322
xAOD::TruthVertex_v1::z
float z() const
Vertex longitudinal distance along the beam line form the origin.
xAOD::TruthParticle_v1::setPz
void setPz(float value)
Set the z component of the particle's momentum.
xAOD::TruthParticle_v1::setDecayVtxLink
void setDecayVtxLink(const ElementLink< TruthVertexContainer > &link)
Set the decay vertex of the particle.
ConstAccessor.h
Helper class to provide constant type-safe access to aux data.
HepMC::status
int status(const T &p)
Definition: MagicNumbers.h:130
TruthEventContainer.h
xAOD::TruthParticle_v1::pdgId
int pdgId() const
PDG ID code.
DataVector::size
size_type size() const noexcept
Returns the number of elements in the collection.
xAOD::TruthVertex_v1::outgoingParticle
const TruthParticle_v1 * outgoingParticle(size_t index) const
Get one of the outgoing particles.
Definition: TruthVertex_v1.cxx:121
DerivationFramework::CollectionMakerHelpers::setupTruthParticle
static xAOD::TruthParticle * setupTruthParticle(const xAOD::TruthParticle &oldPart, xAOD::TruthParticleContainer *partCont)
Definition: CollectionMakerHelpers.cxx:84
xAOD::TruthParticle_v1::m
virtual double m() const override final
The mass of the particle.
HepMCHelpers.h
DerivationFramework::CollectionMakerHelpers::addTruthVertex
static int addTruthVertex(const xAOD::TruthVertex &oldVert, xAOD::TruthParticleContainer *partCont, xAOD::TruthVertexContainer *vertCont, std::vector< int > &seenParticles, const int generations)
Helper function to set up a truth particle based on an old particle.
Definition: CollectionMakerHelpers.cxx:25
xAOD::TruthVertex_v1::setX
void setX(float value)
Set the x displacement of the vertex.
SG::ConstAccessor::withDefault
const_reference_type withDefault(const ELT &e, const T &deflt) const
Fetch the variable for one element, as a const reference, with a default.