ATLAS Offline Software
MultiLeptonFilter.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3 */
4 
7 
8 
9 MultiLeptonFilter::MultiLeptonFilter(const std::string& name, ISvcLocator* pSvcLocator)
10 : GenFilter(name, pSvcLocator)
11 {
12  declareProperty("Ptcut", m_Ptmin = 10000.);
13  declareProperty("Etacut", m_EtaRange = 10.0);
14  declareProperty("NLeptons", m_NLeptons = 4);
15 }
16 
17 
20  int numLeptons = 0;
21  for (itr = events()->begin(); itr != events()->end(); ++itr) {
22  const HepMC::GenEvent* genEvt = *itr;
23  for (const auto& part: *genEvt) {
24  if ( !MC::isStable(part)) continue;
25  if ( MC::isElectron(part) || MC::isMuon(part) ) {
26  if (part->momentum().perp() >= m_Ptmin && std::abs(part->momentum().pseudoRapidity()) <= m_EtaRange) {
27  numLeptons += 1;
28  }
29  }
30  }
31  }
32  ATH_MSG_DEBUG("Found " << numLeptons << " Leptons");
33  setFilterPassed(numLeptons >= m_NLeptons);
34  return StatusCode::SUCCESS;
35 }
LArG4FSStartPointFilter.part
part
Definition: LArG4FSStartPointFilter.py:21
DataModel_detail::const_iterator
Const iterator class for DataVector/DataList.
Definition: DVLIterator.h:82
MultiLeptonFilter::m_Ptmin
double m_Ptmin
Definition: MultiLeptonFilter.h:18
AthCommonDataStore< AthCommonMsg< Algorithm > >::declareProperty
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T > &t)
Definition: AthCommonDataStore.h:145
PlotCalibFromCool.begin
begin
Definition: PlotCalibFromCool.py:94
python.DataFormatRates.events
events
Definition: DataFormatRates.py:105
MultiLeptonFilter::MultiLeptonFilter
MultiLeptonFilter(const std::string &name, ISvcLocator *pSvcLocator)
Definition: MultiLeptonFilter.cxx:9
GenFilter
Base class for event generator filtering modules.
Definition: GenFilter.h:30
MultiLeptonFilter.h
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::EgammaHelpers::isElectron
bool isElectron(const xAOD::Egamma *eg)
is the object an electron (not Fwd)
Definition: EgammaxAODHelpers.cxx:12
MultiLeptonFilter::m_NLeptons
int m_NLeptons
Definition: MultiLeptonFilter.h:20
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:195
MC::isStable
bool isStable(const T &p)
Definition: HepMCHelpers.h:30
MultiLeptonFilter::m_EtaRange
double m_EtaRange
Definition: MultiLeptonFilter.h:19
HepMCHelpers.h
MultiLeptonFilter::filterEvent
virtual StatusCode filterEvent()
Definition: MultiLeptonFilter.cxx:18
isMuon
bool isMuon(const T &p)
Definition: AtlasPID.h:145