ATLAS Offline Software
xAODLeptonFilter.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2021 CERN for the benefit of the ATLAS collaboration
3 */
4 
6 #include <cmath>
8 
9 xAODLeptonFilter::xAODLeptonFilter(const std::string& name, ISvcLocator* pSvcLocator)
10  : GenFilter(name,pSvcLocator)
11 {
12 
13 }
14 
15 
17 
18  // Retrieve TruthElectrons container
19  const xAOD::TruthParticleContainer* xTruthParticleContainerElectron;
20  if (evtStore()->retrieve(xTruthParticleContainerElectron, "TruthElectrons").isFailure()) {
21  ATH_MSG_ERROR("No TruthParticle collection with name " << "TruthElectrons" << " found in StoreGate!");
22  return StatusCode::FAILURE;
23  }
24  // Retrieve TruthMuons container
25  const xAOD::TruthParticleContainer* xTruthParticleContainerMuon;
26  if (evtStore()->retrieve(xTruthParticleContainerMuon, "TruthMuons").isFailure()) {
27  ATH_MSG_ERROR("No TruthParticle collection with name " << "TruthMuons" << " found in StoreGate!");
28  return StatusCode::FAILURE;
29  }
30 
31  double leading_lepton_pt_e = 0;
32  double leading_lepton_pt_mu = 0;
33  double leading_lepton_pt = 0;
34 
35 
36  // Loop over xTruthParticleContainerElectron
37  unsigned int nParticlesElectrons = xTruthParticleContainerElectron->size();
38  for (unsigned int iPart=0; iPart<nParticlesElectrons; ++iPart) {
39  const xAOD::TruthParticle* part = (*xTruthParticleContainerElectron)[iPart];
40 
41 
42 
43 
44  if (MC::isStable(part) && MC::isElectron(part)){ //electron
45  const double pT = part->pt();
46  const double eta = part->abseta();
47  if (pT > leading_lepton_pt_e && std::abs(eta) <= m_EtaRange) {
48  leading_lepton_pt_e = pT;
49  }
50  }
51  }
52 
53  // Loop over xTruthParticleContainerMuon
54  unsigned int nParticlesMuons = xTruthParticleContainerMuon->size();
55  for (unsigned int iPart=0; iPart<nParticlesMuons; ++iPart) {
56  const xAOD::TruthParticle* part = (*xTruthParticleContainerMuon)[iPart];
57 
58  if (MC::isStable(part) && MC::isMuon(part)){ //Muon
59  const double pT = part->pt();
60  const double eta = part->abseta();
61  if (pT > leading_lepton_pt_mu && std::abs(eta) <= m_EtaRange) {
62  leading_lepton_pt_mu = pT;
63  }
64  }
65  }
66 
67  if (leading_lepton_pt_e > leading_lepton_pt_mu){
68  leading_lepton_pt = leading_lepton_pt_e;
69  } else {
70  leading_lepton_pt = leading_lepton_pt_mu;
71  }
72 
73 
74  ATH_MSG_DEBUG ( "Leading lepton pt = " << leading_lepton_pt << "within |eta| <= " << m_EtaRange);
75 
76  if (leading_lepton_pt < m_Ptmin) {
77  setFilterPassed(false);
78  ATH_MSG_DEBUG( "Fail: no e or mu found "
79  << " with pT >= " << m_Ptmin);
80  } else if (leading_lepton_pt >= m_Ptmax) {
81  setFilterPassed(false);
82  ATH_MSG_DEBUG ( "Fail: high pt lepton veto "
83  << " pT < " << m_Ptmax );
84  } else {
85  setFilterPassed(true);
86  ATH_MSG_DEBUG ( "Within min and max pt cuts " << m_Ptmin << ", "
87  << m_Ptmax );
88  }
89  return StatusCode::SUCCESS;
90 }
LArG4FSStartPointFilter.part
part
Definition: LArG4FSStartPointFilter.py:21
python.PyKernel.retrieve
def retrieve(aClass, aKey=None)
Definition: PyKernel.py:110
CalculateHighPtTerm.pT
pT
Definition: ICHEP2016/CalculateHighPtTerm.py:57
eta
Scalar eta() const
pseudorapidity method
Definition: AmgMatrixBasePlugin.h:79
xAODLeptonFilter::m_EtaRange
Gaudi::Property< double > m_EtaRange
Definition: xAODLeptonFilter.h:52
xAODLeptonFilter::xAODLeptonFilter
xAODLeptonFilter(const std::string &name, ISvcLocator *pSvcLocator)
Constructor.
Definition: xAODLeptonFilter.cxx:9
xAODLeptonFilter::m_Ptmax
Gaudi::Property< double > m_Ptmax
Definition: xAODLeptonFilter.h:51
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
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
xAODLeptonFilter::filterEvent
virtual StatusCode filterEvent()
Do the filtering.
Definition: xAODLeptonFilter.cxx:16
xAOD::EgammaHelpers::isElectron
bool isElectron(const xAOD::Egamma *eg)
is the object an electron (not Fwd)
Definition: EgammaxAODHelpers.cxx:12
xAOD::TruthParticle_v1
Class describing a truth particle in the MC record.
Definition: TruthParticle_v1.h:41
DataVector
Derived DataVector<T>.
Definition: DataVector.h:581
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:195
xAODLeptonFilter.h
MC::isStable
bool isStable(const T &p)
Definition: HepMCHelpers.h:30
xAODLeptonFilter::m_Ptmin
Gaudi::Property< double > m_Ptmin
Definition: xAODLeptonFilter.h:50
DataVector::size
size_type size() const noexcept
Returns the number of elements in the collection.
HepMCHelpers.h
isMuon
bool isMuon(const T &p)
Definition: AtlasPID.h:145