ATLAS Offline Software
xAODSameParticleHardScatteringFilter.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 
8  : GenFilter(name, pSvcLocator)
9 {
10  declareProperty("PDGParent", m_PDGParent);
11  declareProperty("PDGChild", m_PDGChild);
12 }
13 
15 {
16  if (m_PDGParent.size() == 0)
17  ATH_MSG_ERROR("PDGParent[] not set ");
18  if (m_PDGChild.size() == 0)
19  ATH_MSG_ERROR("PDGChild[] not set ");
20  for (int i = 0; i < int(m_PDGParent.size()); i++)
21  ATH_MSG_DEBUG("PDGParent[" << i << "] = " << m_PDGParent[i]);
22  for (int i = 0; i < int(m_PDGChild.size()); i++)
23  ATH_MSG_DEBUG("PDGChild[" << i << "] = " << m_PDGChild[i]);
24  return StatusCode::SUCCESS;
25 }
26 
28 {
29  ATH_MSG_DEBUG(" SameParticleHardScattering filtering for: Parent --> " << m_PDGParent[0]
30  << " and parent " << -m_PDGParent[0]
31  << ", Child --> " << m_PDGChild[0]);
32  int N_Parent[2];
33  N_Parent[0] = 0;
34  N_Parent[1] = 0;
35 
36 // Retrieve TruthGen container from xAOD Gen slimmer, contains all particles witout barcode_zero and
37 // duplicated barcode ones
38  const xAOD::TruthParticleContainer* xTruthParticleContainer;
39  if (evtStore()->retrieve(xTruthParticleContainer, "TruthGen").isFailure()) {
40  ATH_MSG_ERROR("No TruthParticle collection with name " << "TruthGen" << " found in StoreGate!");
41  return StatusCode::FAILURE;
42  }
43 
44  // Loop over all particles in the event
45  unsigned int nPart = xTruthParticleContainer->size();
46  for (unsigned int iPart = 0; iPart < nPart; ++iPart) {
47  const xAOD::TruthParticle* part = (*xTruthParticleContainer)[iPart];
48 
49  int id = part->pdgId();
50  if (std::abs(id) != m_PDGChild[0])
51  continue; // searching for only b-quarks
52 
53  // a pointer to the production vertex
54  const xAOD::TruthVertex* productionVtx = part->prodVtx();
55 
56  // Verify if we got a valid pointer and retrieve the number of parents
57  if (!productionVtx)
58  continue;
59  // Incoming particle range check
60 
61  if (productionVtx->nIncomingParticles() < 2) // we are looking for excited tau-leptons produced in b-quark b-antiquark scattering
62  continue;
63  for(size_t thisParent_id=0; thisParent_id < part->prodVtx()->nIncomingParticles(); thisParent_id++)
64  {
65  auto thisParent = part->prodVtx()->incomingParticle(thisParent_id);
66  ATH_MSG_DEBUG(" SelectBQuarkScattering Filter: parent ==> " << thisParent->pdgId() << " child ===> " << part->pdgId());
67  if (thisParent->pdgId() == m_PDGParent[0])
68  {
69  N_Parent[0]++;
70  }
71  if (thisParent->pdgId() == -m_PDGParent[0])
72  {
73  N_Parent[1]++;
74  }
75  }
76  } //loop over TruthParticles
77 
78  setFilterPassed(N_Parent[0] >= 1 && N_Parent[1] >= 1);
79  return StatusCode::SUCCESS;
80 }
LArG4FSStartPointFilter.part
part
Definition: LArG4FSStartPointFilter.py:21
python.PyKernel.retrieve
def retrieve(aClass, aKey=None)
Definition: PyKernel.py:110
xAODSameParticleHardScatteringFilter::m_PDGParent
std::vector< int > m_PDGParent
Definition: xAODSameParticleHardScatteringFilter.h:25
xAODSameParticleHardScatteringFilter::filterInitialize
virtual StatusCode filterInitialize()
Definition: xAODSameParticleHardScatteringFilter.cxx:14
CaloCellPos2Ntuple.int
int
Definition: CaloCellPos2Ntuple.py:24
AthCommonDataStore< AthCommonMsg< Algorithm > >::declareProperty
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T > &t)
Definition: AthCommonDataStore.h:145
xAODSameParticleHardScatteringFilter::filterEvent
virtual StatusCode filterEvent()
Definition: xAODSameParticleHardScatteringFilter.cxx:27
AthCommonDataStore< AthCommonMsg< Algorithm > >::evtStore
ServiceHandle< StoreGateSvc > & evtStore()
The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.
Definition: AthCommonDataStore.h:85
GenFilter
Base class for event generator filtering modules.
Definition: GenFilter.h:30
ATH_MSG_ERROR
#define ATH_MSG_ERROR(x)
Definition: AthMsgStreamMacros.h:33
xAODSameParticleHardScatteringFilter::xAODSameParticleHardScatteringFilter
xAODSameParticleHardScatteringFilter(const std::string &name, ISvcLocator *pSvcLocator)
Definition: xAODSameParticleHardScatteringFilter.cxx:7
lumiFormat.i
int i
Definition: lumiFormat.py:92
EL::StatusCode
::StatusCode StatusCode
StatusCode definition for legacy code.
Definition: PhysicsAnalysis/D3PDTools/EventLoop/EventLoop/StatusCode.h:22
ATH_MSG_DEBUG
#define ATH_MSG_DEBUG(x)
Definition: AthMsgStreamMacros.h:29
xAOD::TruthParticle_v1
Class describing a truth particle in the MC record.
Definition: TruthParticle_v1.h:41
xAODSameParticleHardScatteringFilter::m_PDGChild
std::vector< int > m_PDGChild
Definition: xAODSameParticleHardScatteringFilter.h:26
DataVector
Derived DataVector<T>.
Definition: DataVector.h:581
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
xAOD::TruthVertex_v1::nIncomingParticles
size_t nIncomingParticles() const
Get the number of incoming particles.
Definition: TruthVertex_v1.cxx:49
xAODSameParticleHardScatteringFilter.h
DataVector::size
size_type size() const noexcept
Returns the number of elements in the collection.