ATLAS Offline Software
xAODParticleFilter.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 // Allows the user to search for particles with specified kinematics.
6 // It will pass if there is an particle with pT and eta or E in the specified
7 // range
9 
10 #include "xAODTruth/TruthEvent.h"
12 
14  ISvcLocator *pSvcLocator)
15  : GenFilter(name, pSvcLocator) {}
16 
18 {
19  ATH_MSG_INFO("Ptcut=" << m_Ptmin);
20  ATH_MSG_INFO("Etacut=" << m_EtaRange);
21  ATH_MSG_INFO("Energycut=" << m_EnergyRange);
22  ATH_MSG_INFO("PDG=" << m_PDGID);
23  ATH_MSG_INFO("MinParts=" << m_MinParts);
24  ATH_MSG_INFO("Exclusive=" << m_Exclusive);
25  return StatusCode::SUCCESS;
26 }
27 
29 {
30  int nParts = 0;
31 
32 // Retrieve TruthGen container from xAOD Gen slimmer, contains all particles witout barcode_zero and
33 // duplicated barcode ones
34  const xAOD::TruthParticleContainer* xTruthParticleContainer;
35  if (evtStore()->retrieve(xTruthParticleContainer, "TruthGen").isFailure()) {
36  ATH_MSG_ERROR("No TruthParticle collection with name " << "TruthGen" << " found in StoreGate!");
37  return StatusCode::FAILURE;
38  }
39 
40  // Loop over all particles in the event and build up the grid
41  unsigned int nPart = xTruthParticleContainer->size();
42  for (unsigned int iPart = 0; iPart < nPart; ++iPart) {
43  const xAOD::TruthParticle* pitr = (*xTruthParticleContainer)[iPart];
44 
45  if (std::abs(pitr->pdgId()) != m_PDGID)
46  continue;
47  if (pitr->pt() >= m_Ptmin && std::abs(pitr->eta()) <= m_EtaRange &&
48  pitr->e() <= m_EnergyRange)
49  {
50  if ((!m_Exclusive) && (m_MinParts == 1))
51  return StatusCode::SUCCESS; // Found at least one particle and we have an inclusive requirement
52 
53  // Count only particles not decaying to themselves
54  bool notSelfDecay = true;
55  if (pitr->decayVtx())
56  {
57  const xAOD::TruthVertex *decayVertex = pitr->decayVtx();
58  int outgoing_particles = decayVertex->nOutgoingParticles();
59  for (int part = 0; part < outgoing_particles;
60  part++)
61  {
62  const xAOD::TruthParticle *child =
63  decayVertex->outgoingParticle(part);
64  if (child->pdgId() != pitr->pdgId())
65  continue;
66  if (child == pitr)
67  continue;
69  continue;
70  notSelfDecay = false;
71  break;
72  }
73  }
74  if (notSelfDecay)
75  nParts++;
76  }
77  }
78  if (m_Exclusive)
79  {
80  setFilterPassed(nParts == m_MinParts);
81  }
82  else
83  {
84  setFilterPassed(nParts >= m_MinParts);
85  }
86  return StatusCode::SUCCESS;
87 }
LArG4FSStartPointFilter.part
part
Definition: LArG4FSStartPointFilter.py:21
python.PyKernel.retrieve
def retrieve(aClass, aKey=None)
Definition: PyKernel.py:110
xAOD::TruthVertex_v1::nOutgoingParticles
size_t nOutgoingParticles() const
Get the number of outgoing particles.
xAODParticleFilter.h
ATH_MSG_INFO
#define ATH_MSG_INFO(x)
Definition: AthMsgStreamMacros.h:31
xAODParticleFilter::m_PDGID
Gaudi::Property< int > m_PDGID
Definition: xAODParticleFilter.h:26
xAODParticleFilter::filterEvent
virtual StatusCode filterEvent() override
Definition: xAODParticleFilter.cxx:28
xAODParticleFilter::filterInitialize
virtual StatusCode filterInitialize() override
Definition: xAODParticleFilter.cxx:17
AthCommonDataStore< AthCommonMsg< Algorithm > >::evtStore
ServiceHandle< StoreGateSvc > & evtStore()
The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.
Definition: AthCommonDataStore.h:85
xAODParticleFilter::m_MinParts
Gaudi::Property< int > m_MinParts
Definition: xAODParticleFilter.h:27
GenFilter
Base class for event generator filtering modules.
Definition: GenFilter.h:30
ATH_MSG_ERROR
#define ATH_MSG_ERROR(x)
Definition: AthMsgStreamMacros.h:33
xAOD::TruthParticle_v1::e
virtual double e() const override final
The total energy of the particle.
HepMC::is_simulation_particle
bool is_simulation_particle(const T &p)
Method to establish if a particle (or barcode) was created during the simulation (TODO update to be s...
Definition: MagicNumbers.h:299
EL::StatusCode
::StatusCode StatusCode
StatusCode definition for legacy code.
Definition: PhysicsAnalysis/D3PDTools/EventLoop/EventLoop/StatusCode.h:22
xAOD::TruthParticle_v1
Class describing a truth particle in the MC record.
Definition: TruthParticle_v1.h:41
xAODParticleFilter::m_Exclusive
Gaudi::Property< bool > m_Exclusive
Definition: xAODParticleFilter.h:28
DataVector
Derived DataVector<T>.
Definition: DataVector.h:581
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
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:195
xAODParticleFilter::xAODParticleFilter
xAODParticleFilter(const std::string &name, ISvcLocator *pSvcLocator)
Definition: xAODParticleFilter.cxx:13
xAOD::TruthParticle_v1::eta
virtual double eta() const override final
The pseudorapidity ( ) of the particle.
Definition: TruthParticle_v1.cxx:174
xAOD::TruthParticle_v1::pt
virtual double pt() const override final
The transverse momentum ( ) of the particle.
Definition: TruthParticle_v1.cxx:166
xAODParticleFilter::m_EtaRange
Gaudi::Property< double > m_EtaRange
Definition: xAODParticleFilter.h:24
xAODParticleFilter::m_EnergyRange
Gaudi::Property< double > m_EnergyRange
Definition: xAODParticleFilter.h:25
TruthEventContainer.h
xAOD::TruthParticle_v1::pdgId
int pdgId() const
PDG ID code.
xAODParticleFilter::m_Ptmin
Gaudi::Property< double > m_Ptmin
Definition: xAODParticleFilter.h:23
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
TruthEvent.h