#include <Configuration/Skimming/src/LeptonRecoSkim.cc>
Description: [one line class summary]
Implementation: [Notes on implementation]
Definition at line 59 of file LeptonRecoSkim.h.
LeptonRecoSkim::LeptonRecoSkim | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 32 of file LeptonRecoSkim.cc.
References NeventsFiltered, NeventsTotal, NHltDiMu3, NHltMu9, NmvaElectrons, and NtotalElectrons.
: hltLabel(iConfig.getParameter<edm::InputTag>("HltLabel")), filterName(iConfig.getParameter<std::string>("@module_label")), m_electronSrc(iConfig.getParameter<edm::InputTag>("electronCollection")), m_pfelectronSrc(iConfig.getParameter<edm::InputTag>("pfElectronCollection")), m_muonSrc(iConfig.getParameter<edm::InputTag>("muonCollection")), m_jetsSrc(iConfig.getParameter<edm::InputTag>("caloJetCollection")), m_pfjetsSrc(iConfig.getParameter<edm::InputTag>("PFJetCollection")), m_ebRecHitsSrc(iConfig.getParameter<edm::InputTag>("ecalBarrelRecHitsCollection")), m_eeRecHitsSrc(iConfig.getParameter<edm::InputTag>("ecalEndcapRecHitsCollection")), useElectronSelection(iConfig.getParameter<bool>("UseElectronSelection")), usePfElectronSelection(iConfig.getParameter<bool>("UsePfElectronSelection")), useMuonSelection(iConfig.getParameter<bool>("UseMuonSelection")), useHtSelection(iConfig.getParameter<bool>("UseHtSelection")), usePFHtSelection(iConfig.getParameter<bool>("UsePFHtSelection")), ptElecMin(iConfig.getParameter<double>("electronPtMin")), ptPfElecMin(iConfig.getParameter<double>("pfElectronPtMin")), nSelectedElectrons(iConfig.getParameter<int>("electronN")), nSelectedPfElectrons(iConfig.getParameter<int>("pfElectronN")), ptGlobalMuonMin(iConfig.getParameter<double>("globalMuonPtMin")), ptTrackerMuonMin(iConfig.getParameter<double>("trackerMuonPtMin")), nSelectedMuons(iConfig.getParameter<int>("muonN")), htMin(iConfig.getParameter<double>("HtMin")), pfHtMin(iConfig.getParameter<double>("PFHtMin")), htJetThreshold(iConfig.getParameter<double>("HtJetThreshold")), pfHtJetThreshold(iConfig.getParameter<double>("PFHtJetThreshold")) { NeventsTotal = 0; NeventsFiltered = 0; NHltMu9 = 0; NHltDiMu3 = 0; NtotalElectrons = 0; NmvaElectrons = 0; }
LeptonRecoSkim::~LeptonRecoSkim | ( | ) |
Definition at line 69 of file LeptonRecoSkim.cc.
{ // do anything here that needs to be done at desctruction time // (e.g. close files, deallocate resources etc.) }
void LeptonRecoSkim::beginJob | ( | void | ) | [private, virtual] |
Reimplemented from edm::EDFilter.
Definition at line 214 of file LeptonRecoSkim.cc.
References firstEvent.
{ firstEvent = true; }
void LeptonRecoSkim::endJob | ( | void | ) | [private, virtual] |
Reimplemented from edm::EDFilter.
Definition at line 221 of file LeptonRecoSkim.cc.
References gather_cfg::cout, filterName, NeventsFiltered, NeventsTotal, NHltDiMu3, NHltMu9, NmvaElectrons, and NtotalElectrons.
{ cout << "Filter Name = " << filterName << endl; cout << "Total number of events = " << NeventsTotal << endl; cout << "Total HLT_Mu9 = " << NHltMu9 << endl; cout << "Total HLT_DoubleMu3 = " << NHltDiMu3 << endl; cout << "Filtered events = " << NeventsFiltered << endl; cout << "Filter Efficiency = " << (float)NeventsFiltered / (float)NeventsTotal << endl; cout << endl; cout << "N total electrons = " << NtotalElectrons << endl; cout << "N mva>-0.1 electrons = " << NmvaElectrons << endl; cout << endl; cout << endl; }
bool LeptonRecoSkim::filter | ( | edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [private, virtual] |
Implements edm::EDFilter.
Definition at line 84 of file LeptonRecoSkim.cc.
References accept(), reco::PFCandidate::e, metsig::electron, ElectronCutPassed, firstEvent, reco::Muon::globalTrack(), reco::PFCandidate::gsfTrackRef(), handleObjects(), HtCutPassed, htJetThreshold, htMin, i, reco::Muon::innerTrack(), reco::Muon::isGlobalMuon(), edm::Ref< C, T, F >::isNull(), reco::Muon::isTrackerMuon(), LogDebug, max(), metsig::muon, MuonCutPassed, reco::GsfElectron::mva(), NeventsFiltered, NeventsTotal, NmvaElectrons, nSelectedElectrons, nSelectedMuons, nSelectedPfElectrons, NtotalElectrons, reco::PFCandidate::particleId(), PfElectronCutPassed, PFHtCutPassed, pfHtJetThreshold, pfHtMin, reco::LeafCandidate::pt(), ptElecMin, ptGlobalMuonMin, ptPfElecMin, ptTrackerMuonMin, theCaloJetCollection, theElectronCollection, theMuonCollection, thePfCandidateCollection, thePFJetCollection, useElectronSelection, useHtSelection, useMuonSelection, usePfElectronSelection, and usePFHtSelection.
{ bool accept = false; NeventsTotal++; ElectronCutPassed = false; PfElectronCutPassed = false; MuonCutPassed = false; HtCutPassed = false; PFHtCutPassed = false; // edm::Handle<TriggerResults> trhv; // iEvent.getByLabel(hltLabel,trhv); // const edm::TriggerNames& triggerNames_ = iEvent.triggerNames(*trhv); // if(trhv->at(triggerNames_.triggerIndex("HLT_Mu9")).accept()) NHltMu9++; // if(trhv->at(triggerNames_.triggerIndex("HLT_DoubleMu3")).accept()) NHltDiMu3++; this->handleObjects(iEvent,iSetup); if(useElectronSelection) { int nElecPassingCut = 0; for(unsigned int i=0; i<theElectronCollection->size(); i++) { GsfElectron electron = (*theElectronCollection)[i]; // if (electron.ecalDrivenSeed()) { float elpt = electron.pt(); // if (electron.sigmaIetaIeta() < 0.002) continue; if(elpt>ptElecMin) { NtotalElectrons++; nElecPassingCut++; if(electron.mva()>-0.1) NmvaElectrons++; } LogDebug("LeptonRecoSkim") << "elpt = " << elpt << endl; // } // closes if (electron.ecalDrivenSeed()) { } if(nElecPassingCut>=nSelectedElectrons) ElectronCutPassed = true; } else ElectronCutPassed = true; if(usePfElectronSelection) { int nPfElecPassingCut = 0; for(unsigned int i=0; i<thePfCandidateCollection->size(); i++) { const reco::PFCandidate& thePfCandidate = (*thePfCandidateCollection)[i]; if(thePfCandidate.particleId()!=reco::PFCandidate::e) continue; if(thePfCandidate.gsfTrackRef().isNull()) continue; float pfelpt = thePfCandidate.pt(); // if (electron.sigmaIetaIeta() < 0.002) continue; if(pfelpt>ptPfElecMin) nPfElecPassingCut++; LogDebug("LeptonRecoSkim") << "pfelpt = " << pfelpt << endl; } if(nPfElecPassingCut>=nSelectedPfElectrons) PfElectronCutPassed = true; } else PfElectronCutPassed = true; if(useMuonSelection) { int nMuonPassingCut = 0; for(unsigned int i=0; i<theMuonCollection->size(); i++) { Muon muon = (*theMuonCollection)[i]; if(! (muon.isGlobalMuon() || muon.isTrackerMuon()) ) continue; const TrackRef siTrack = muon.innerTrack(); const TrackRef globalTrack = muon.globalTrack(); float muonpt; float ptMuonMin; // if (siTrack.isNonnull() && globalTrack.isNonnull()) { if (muon.isGlobalMuon() && muon.isTrackerMuon()) { muonpt = max(siTrack->pt(), globalTrack->pt()); ptMuonMin = ptGlobalMuonMin; } else if (muon.isGlobalMuon() && !(muon.isTrackerMuon())) { muonpt = globalTrack->pt(); ptMuonMin = ptGlobalMuonMin; } else if (muon.isTrackerMuon() && !(muon.isGlobalMuon())) { muonpt = siTrack->pt(); ptMuonMin = ptTrackerMuonMin; } else { muonpt = 0; // if we get here this is a STA only muon ptMuonMin = 999999.; } if(muonpt>ptMuonMin) nMuonPassingCut++; LogDebug("RecoSelectorCuts") << "muonpt = " << muonpt << endl; } if(nMuonPassingCut>=nSelectedMuons) MuonCutPassed = true; } else MuonCutPassed = true; if(useHtSelection) { double Ht = 0; for(unsigned int i=0; i<theCaloJetCollection->size(); i++) { // if((*theCaloJetCollection)[i].eta()<2.6 && (*theCaloJetCollection)[i].emEnergyFraction() <= 0.01) continue; if((*theCaloJetCollection)[i].pt()>htJetThreshold) Ht += (*theCaloJetCollection)[i].pt(); } if(Ht>htMin) HtCutPassed = true; } else HtCutPassed = true; if(usePFHtSelection) { double PFHt = 0; for(unsigned int i=0; i<thePFJetCollection->size(); i++) { if((*thePFJetCollection)[i].pt()>pfHtJetThreshold) PFHt += (*thePFJetCollection)[i].pt(); } if(PFHt>pfHtMin) PFHtCutPassed = true; } else PFHtCutPassed = true; if(PfElectronCutPassed && ElectronCutPassed && MuonCutPassed && HtCutPassed && PFHtCutPassed ) accept = true; if(accept) NeventsFiltered++; firstEvent = false; return accept; }
void LeptonRecoSkim::handleObjects | ( | const edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [private] |
Definition at line 238 of file LeptonRecoSkim.cc.
References edm::EventSetup::get(), edm::Event::getByLabel(), m_ebRecHitsSrc, m_eeRecHitsSrc, m_electronSrc, m_jetsSrc, m_muonSrc, m_pfelectronSrc, m_pfjetsSrc, edm::ESHandle< T >::product(), edm::Handle< T >::product(), theCaloGeometry, theCaloJetCollection, theCaloTopology, theEcalBarrelCollection, theEcalEndcapCollection, theElectronCollection, theMuonCollection, thePfCandidateCollection, and thePFJetCollection.
Referenced by filter().
{ //Get the electrons Handle<GsfElectronCollection> theElectronCollectionHandle; iEvent.getByLabel(m_electronSrc, theElectronCollectionHandle); theElectronCollection = theElectronCollectionHandle.product(); //Get the pf electrons Handle<reco::PFCandidateCollection> thePfCandidateHandle; iEvent.getByLabel(m_pfelectronSrc, thePfCandidateHandle); thePfCandidateCollection = thePfCandidateHandle.product(); //Get the Muons Handle<MuonCollection> theMuonCollectionHandle; iEvent.getByLabel(m_muonSrc, theMuonCollectionHandle); theMuonCollection = theMuonCollectionHandle.product(); //Get the CaloJets Handle<CaloJetCollection> theCaloJetCollectionHandle; iEvent.getByLabel(m_jetsSrc, theCaloJetCollectionHandle); theCaloJetCollection = theCaloJetCollectionHandle.product(); //Get the PfJets Handle<PFJetCollection> thePFJetCollectionHandle; iEvent.getByLabel(m_pfjetsSrc, thePFJetCollectionHandle); thePFJetCollection = thePFJetCollectionHandle.product(); //Get the ECAL rechhits to clean the spikes // Get EB RecHits edm::Handle<EcalRecHitCollection> ebRecHitsHandle; iEvent.getByLabel(m_ebRecHitsSrc, ebRecHitsHandle); theEcalBarrelCollection = ebRecHitsHandle.product(); // Get EE RecHits edm::Handle<EcalRecHitCollection> eeRecHitsHandle; iEvent.getByLabel(m_eeRecHitsSrc, eeRecHitsHandle); theEcalEndcapCollection = eeRecHitsHandle.product(); // Get topology for spike cleaning edm::ESHandle<CaloGeometry> geometryHandle; iSetup.get<CaloGeometryRecord>().get(geometryHandle); theCaloGeometry = geometryHandle.product(); // theCaloBarrelSubdetTopology = new EcalBarrelTopology(geometryHandle); // theCaloEndcapSubdetTopology = new EcalEndcapTopology(geometryHandle); edm::ESHandle<CaloTopology> pTopology; iSetup.get<CaloTopologyRecord>().get(pTopology); theCaloTopology = pTopology.product(); }
bool LeptonRecoSkim::ElectronCutPassed [private] |
Definition at line 104 of file LeptonRecoSkim.h.
Referenced by filter().
std::string LeptonRecoSkim::filterName [private] |
Definition at line 74 of file LeptonRecoSkim.h.
Referenced by endJob().
bool LeptonRecoSkim::firstEvent [private] |
Definition at line 93 of file LeptonRecoSkim.h.
Referenced by beginJob(), and filter().
edm::InputTag LeptonRecoSkim::hltLabel [private] |
Definition at line 73 of file LeptonRecoSkim.h.
bool LeptonRecoSkim::HtCutPassed [private] |
Definition at line 107 of file LeptonRecoSkim.h.
Referenced by filter().
double LeptonRecoSkim::htJetThreshold [private] |
Definition at line 120 of file LeptonRecoSkim.h.
Referenced by filter().
double LeptonRecoSkim::htMin [private] |
Definition at line 118 of file LeptonRecoSkim.h.
Referenced by filter().
edm::InputTag LeptonRecoSkim::m_ebRecHitsSrc [private] |
Definition at line 80 of file LeptonRecoSkim.h.
Referenced by handleObjects().
edm::InputTag LeptonRecoSkim::m_eeRecHitsSrc [private] |
Definition at line 81 of file LeptonRecoSkim.h.
Referenced by handleObjects().
edm::InputTag LeptonRecoSkim::m_electronSrc [private] |
Definition at line 75 of file LeptonRecoSkim.h.
Referenced by handleObjects().
edm::InputTag LeptonRecoSkim::m_jetsSrc [private] |
Definition at line 78 of file LeptonRecoSkim.h.
Referenced by handleObjects().
edm::InputTag LeptonRecoSkim::m_muonSrc [private] |
Definition at line 77 of file LeptonRecoSkim.h.
Referenced by handleObjects().
edm::InputTag LeptonRecoSkim::m_pfelectronSrc [private] |
Definition at line 76 of file LeptonRecoSkim.h.
Referenced by handleObjects().
edm::InputTag LeptonRecoSkim::m_pfjetsSrc [private] |
Definition at line 79 of file LeptonRecoSkim.h.
Referenced by handleObjects().
bool LeptonRecoSkim::MuonCutPassed [private] |
Definition at line 106 of file LeptonRecoSkim.h.
Referenced by filter().
int LeptonRecoSkim::NeventsFiltered [private] |
Definition at line 126 of file LeptonRecoSkim.h.
Referenced by endJob(), filter(), and LeptonRecoSkim().
int LeptonRecoSkim::NeventsTotal [private] |
Definition at line 125 of file LeptonRecoSkim.h.
Referenced by endJob(), filter(), and LeptonRecoSkim().
int LeptonRecoSkim::NHltDiMu3 [private] |
Definition at line 128 of file LeptonRecoSkim.h.
Referenced by endJob(), and LeptonRecoSkim().
int LeptonRecoSkim::NHltMu9 [private] |
Definition at line 127 of file LeptonRecoSkim.h.
Referenced by endJob(), and LeptonRecoSkim().
int LeptonRecoSkim::NmvaElectrons [private] |
Definition at line 131 of file LeptonRecoSkim.h.
Referenced by endJob(), filter(), and LeptonRecoSkim().
int LeptonRecoSkim::nSelectedElectrons [private] |
Definition at line 113 of file LeptonRecoSkim.h.
Referenced by filter().
int LeptonRecoSkim::nSelectedMuons [private] |
Definition at line 117 of file LeptonRecoSkim.h.
Referenced by filter().
int LeptonRecoSkim::nSelectedPfElectrons [private] |
Definition at line 114 of file LeptonRecoSkim.h.
Referenced by filter().
int LeptonRecoSkim::NtotalElectrons [private] |
Definition at line 130 of file LeptonRecoSkim.h.
Referenced by endJob(), filter(), and LeptonRecoSkim().
bool LeptonRecoSkim::PfElectronCutPassed [private] |
Definition at line 105 of file LeptonRecoSkim.h.
Referenced by filter().
bool LeptonRecoSkim::PFHtCutPassed [private] |
Definition at line 108 of file LeptonRecoSkim.h.
Referenced by filter().
double LeptonRecoSkim::pfHtJetThreshold [private] |
Definition at line 121 of file LeptonRecoSkim.h.
Referenced by filter().
double LeptonRecoSkim::pfHtMin [private] |
Definition at line 119 of file LeptonRecoSkim.h.
Referenced by filter().
double LeptonRecoSkim::ptElecMin [private] |
Definition at line 111 of file LeptonRecoSkim.h.
Referenced by filter().
double LeptonRecoSkim::ptGlobalMuonMin [private] |
Definition at line 115 of file LeptonRecoSkim.h.
Referenced by filter().
double LeptonRecoSkim::ptPfElecMin [private] |
Definition at line 112 of file LeptonRecoSkim.h.
Referenced by filter().
double LeptonRecoSkim::ptTrackerMuonMin [private] |
Definition at line 116 of file LeptonRecoSkim.h.
Referenced by filter().
const CaloGeometry* LeptonRecoSkim::theCaloGeometry [private] |
Definition at line 91 of file LeptonRecoSkim.h.
Referenced by handleObjects().
const reco::CaloJetCollection* LeptonRecoSkim::theCaloJetCollection [private] |
Definition at line 86 of file LeptonRecoSkim.h.
Referenced by filter(), and handleObjects().
const CaloTopology* LeptonRecoSkim::theCaloTopology [private] |
Definition at line 90 of file LeptonRecoSkim.h.
Referenced by handleObjects().
const EcalRecHitCollection* LeptonRecoSkim::theEcalBarrelCollection [private] |
Definition at line 88 of file LeptonRecoSkim.h.
Referenced by handleObjects().
const EcalRecHitCollection* LeptonRecoSkim::theEcalEndcapCollection [private] |
Definition at line 89 of file LeptonRecoSkim.h.
Referenced by handleObjects().
const reco::GsfElectronCollection* LeptonRecoSkim::theElectronCollection [private] |
Definition at line 83 of file LeptonRecoSkim.h.
Referenced by filter(), and handleObjects().
const reco::MuonCollection* LeptonRecoSkim::theMuonCollection [private] |
Definition at line 85 of file LeptonRecoSkim.h.
Referenced by filter(), and handleObjects().
const reco::PFCandidateCollection* LeptonRecoSkim::thePfCandidateCollection [private] |
Definition at line 84 of file LeptonRecoSkim.h.
Referenced by filter(), and handleObjects().
const reco::PFJetCollection* LeptonRecoSkim::thePFJetCollection [private] |
Definition at line 87 of file LeptonRecoSkim.h.
Referenced by filter(), and handleObjects().
bool LeptonRecoSkim::useElectronSelection [private] |
Definition at line 97 of file LeptonRecoSkim.h.
Referenced by filter().
bool LeptonRecoSkim::useHtSelection [private] |
Definition at line 100 of file LeptonRecoSkim.h.
Referenced by filter().
bool LeptonRecoSkim::useMuonSelection [private] |
Definition at line 99 of file LeptonRecoSkim.h.
Referenced by filter().
bool LeptonRecoSkim::usePfElectronSelection [private] |
Definition at line 98 of file LeptonRecoSkim.h.
Referenced by filter().
bool LeptonRecoSkim::usePFHtSelection [private] |
Definition at line 101 of file LeptonRecoSkim.h.
Referenced by filter().