ATLAS Offline Software
xAODParentChildFilter.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 
7 xAODParentChildFilter::xAODParentChildFilter(const std::string& name, ISvcLocator* pSvcLocator)
8  : GenFilter(name,pSvcLocator)
9 {
10 
11 }
12 
13 
15  if (m_PDGParent.size() == 0) ATH_MSG_ERROR("PDGParent[] not set ");
16  if (m_PDGChild.size() == 0) ATH_MSG_ERROR("PDGChild[] not set ");
17  for (int i=0; i < int(m_PDGParent.size()); i++) ATH_MSG_INFO("PDGParent["<<i<<"] = " << m_PDGParent[i]);
18  ATH_MSG_INFO("PtMinParent = " << m_PtMinParent);
19  ATH_MSG_INFO("PtMaxParent = " << m_PtMaxParent);
20  ATH_MSG_INFO("MassMinParent = " << m_MassMinParent);
21  ATH_MSG_INFO("MassMaxParent = " << m_MassMaxParent);
22  ATH_MSG_INFO("EtaRangeParent = " << m_EtaRangeParent);
23  for (int i=0; i < int(m_PDGChild.size()); i++) ATH_MSG_INFO("PDGChild["<<i<<"] = " << m_PDGChild[i]);
24  ATH_MSG_INFO("PtMinChild = " << m_PtMinChild);
25  ATH_MSG_INFO("EtaRangeChild = " << m_EtaRangeChild);
26  return StatusCode::SUCCESS;
27 }
28 
29 
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  int okPDGParent = 0;
46  for (int i = 0; i < int(m_PDGParent.size()); i++) if (std::abs(pitr->pdgId()) == m_PDGParent[i]) okPDGParent=1;
47  if ( (m_PDGParent[0] == 0) || (okPDGParent
48  && pitr->pt() >= m_PtMinParent
49  && (pitr->pt() < m_PtMaxParent)
50  && (pitr->m() >= m_MassMinParent)
51  && (pitr->m() < m_MassMaxParent)
52  && std::abs(pitr->eta()) < m_EtaRangeParent)) {
53 
54  // Check if has end_vertex (skips initial protons)
55  if (!pitr->decayVtx()) continue;
56 
57  const xAOD::TruthVertex* decayVertex = pitr->decayVtx();
58  int num_outgoing_particles = decayVertex->nOutgoingParticles();
59 
60  for (int iOutPart = 0; iOutPart< num_outgoing_particles; iOutPart++) {
61  const xAOD::TruthParticle* thisChild = decayVertex->outgoingParticle(iOutPart);
62  int okPDGChild = 0;
63  for (int i=0; i < int(m_PDGChild.size()); i++) if (std::abs(thisChild->pdgId()) == m_PDGChild[i]) okPDGChild=1;
64  if ( thisChild->pdgId() != pitr->pdgId()
65  && (m_PDGChild[0] == 0 || okPDGChild)
66  && thisChild->pt() > m_PtMinChild
67  && std::abs(thisChild->eta()) < m_EtaRangeChild ) {
68  return StatusCode::SUCCESS;
69  }
70  }
71  }
72  }// TruthParticles loop
73 
74  setFilterPassed(false);
75  return StatusCode::SUCCESS;
76 }
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.
xAODParentChildFilter::filterInitialize
virtual StatusCode filterInitialize()
Definition: xAODParentChildFilter.cxx:14
ATH_MSG_INFO
#define ATH_MSG_INFO(x)
Definition: AthMsgStreamMacros.h:31
CaloCellPos2Ntuple.int
int
Definition: CaloCellPos2Ntuple.py:24
xAODParentChildFilter::m_PtMinChild
Gaudi::Property< double > m_PtMinChild
Definition: xAODParentChildFilter.h:32
xAODParentChildFilter::m_MassMinParent
Gaudi::Property< double > m_MassMinParent
Definition: xAODParentChildFilter.h: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
xAODParentChildFilter::m_PtMinParent
Gaudi::Property< double > m_PtMinParent
Definition: xAODParentChildFilter.h:25
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
xAOD::TruthParticle_v1
Class describing a truth particle in the MC record.
Definition: TruthParticle_v1.h:41
xAODParentChildFilter::filterEvent
virtual StatusCode filterEvent()
Definition: xAODParentChildFilter.cxx:30
xAODParentChildFilter.h
xAODParentChildFilter::m_PtMaxParent
Gaudi::Property< double > m_PtMaxParent
Definition: xAODParentChildFilter.h:26
DataVector
Derived DataVector<T>.
Definition: DataVector.h:581
xAODParentChildFilter::m_PDGParent
Gaudi::Property< std::vector< int > > m_PDGParent
Definition: xAODParentChildFilter.h:29
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
xAOD::TruthParticle_v1::eta
virtual double eta() const override final
The pseudorapidity ( ) of the particle.
Definition: TruthParticle_v1.cxx:174
xAODParentChildFilter::m_PDGChild
Gaudi::Property< std::vector< int > > m_PDGChild
Definition: xAODParentChildFilter.h:31
xAODParentChildFilter::m_EtaRangeParent
Gaudi::Property< double > m_EtaRangeParent
Definition: xAODParentChildFilter.h:30
xAODParentChildFilter::xAODParentChildFilter
xAODParentChildFilter(const std::string &name, ISvcLocator *pSvcLocator)
Definition: xAODParentChildFilter.cxx:7
xAOD::TruthParticle_v1::pt
virtual double pt() const override final
The transverse momentum ( ) of the particle.
Definition: TruthParticle_v1.cxx:166
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
xAODParentChildFilter::m_EtaRangeChild
Gaudi::Property< double > m_EtaRangeChild
Definition: xAODParentChildFilter.h:33
xAOD::TruthParticle_v1::m
virtual double m() const override final
The mass of the particle.
xAODParentChildFilter::m_MassMaxParent
Gaudi::Property< double > m_MassMaxParent
Definition: xAODParentChildFilter.h:28