Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:03:41

0001 
0002 // -*- C++ -*-
0003 //
0004 // Package:    LeptonRecoSkim
0005 // Class:      LeptonRecoSkim
0006 //
0007 /**\class LeptonRecoSkim LeptonRecoSkim.cc Configuration/Skimming/src/LeptonRecoSkim.cc
0008 
0009  Description: [one line class summary]
0010 
0011  Implementation:
0012      [Notes on implementation]
0013 */
0014 //
0015 // Original Author:  Massimiliano Chiorboli,40 4-A01,+41227671535,
0016 //         Created:  Wed Mar 31 21:49:08 CEST 2010
0017 // $Id: LeptonRecoSkim.cc,v 1.1 2010/11/05 18:37:51 torimoto Exp $
0018 //
0019 //
0020 
0021 #include "Configuration/Skimming/interface/LeptonRecoSkim.h"
0022 
0023 using namespace edm;
0024 using namespace reco;
0025 using namespace std;
0026 
0027 //
0028 // constructors and destructor
0029 //
0030 LeptonRecoSkim::LeptonRecoSkim(const edm::ParameterSet& iConfig)
0031     : m_CaloGeoToken(esConsumes()),
0032       m_CaloTopoToken(esConsumes()),
0033       hltLabel(iConfig.getParameter<edm::InputTag>("HltLabel")),
0034       filterName(iConfig.getParameter<std::string>("@module_label")),
0035       gsfElectronCollectionToken_(consumes(iConfig.getParameter<edm::InputTag>("electronCollection"))),
0036       pfCandidateCollectionToken_(consumes(iConfig.getParameter<edm::InputTag>("pfElectronCollection"))),
0037       muonCollectionToken_(consumes(iConfig.getParameter<edm::InputTag>("muonCollection"))),
0038       caloJetCollectionToken_(consumes(iConfig.getParameter<edm::InputTag>("caloJetCollection"))),
0039       pfJetCollectionToken_(consumes(iConfig.getParameter<edm::InputTag>("PFJetCollection"))),
0040       ebRecHitCollectionToken_(consumes(iConfig.getParameter<edm::InputTag>("ecalBarrelRecHitsCollection"))),
0041       eeRecHitCollectionToken_(consumes(iConfig.getParameter<edm::InputTag>("ecalEndcapRecHitsCollection"))),
0042       useElectronSelection(iConfig.getParameter<bool>("UseElectronSelection")),
0043       usePfElectronSelection(iConfig.getParameter<bool>("UsePfElectronSelection")),
0044       useMuonSelection(iConfig.getParameter<bool>("UseMuonSelection")),
0045       useHtSelection(iConfig.getParameter<bool>("UseHtSelection")),
0046       usePFHtSelection(iConfig.getParameter<bool>("UsePFHtSelection")),
0047       ptElecMin(iConfig.getParameter<double>("electronPtMin")),
0048       ptPfElecMin(iConfig.getParameter<double>("pfElectronPtMin")),
0049       nSelectedElectrons(iConfig.getParameter<int>("electronN")),
0050       nSelectedPfElectrons(iConfig.getParameter<int>("pfElectronN")),
0051       ptGlobalMuonMin(iConfig.getParameter<double>("globalMuonPtMin")),
0052       ptTrackerMuonMin(iConfig.getParameter<double>("trackerMuonPtMin")),
0053       nSelectedMuons(iConfig.getParameter<int>("muonN")),
0054       htMin(iConfig.getParameter<double>("HtMin")),
0055       pfHtMin(iConfig.getParameter<double>("PFHtMin")),
0056       htJetThreshold(iConfig.getParameter<double>("HtJetThreshold")),
0057       pfHtJetThreshold(iConfig.getParameter<double>("PFHtJetThreshold")) {
0058   NeventsTotal = 0;
0059   NeventsFiltered = 0;
0060   NHltMu9 = 0;
0061   NHltDiMu3 = 0;
0062 
0063   NtotalElectrons = 0;
0064   NmvaElectrons = 0;
0065 }
0066 
0067 LeptonRecoSkim::~LeptonRecoSkim() {
0068   // do anything here that needs to be done at desctruction time
0069   // (e.g. close files, deallocate resources etc.)
0070 }
0071 
0072 //
0073 // member functions
0074 //
0075 
0076 // ------------ method called on each new Event  ------------
0077 bool LeptonRecoSkim::filter(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0078   bool accept = false;
0079 
0080   NeventsTotal++;
0081 
0082   ElectronCutPassed = false;
0083   PfElectronCutPassed = false;
0084   MuonCutPassed = false;
0085   HtCutPassed = false;
0086   PFHtCutPassed = false;
0087 
0088   //  edm::Handle<TriggerResults> trhv;
0089   //  iEvent.getByLabel(hltLabel,trhv);
0090 
0091   //  const edm::TriggerNames& triggerNames_ = iEvent.triggerNames(*trhv);
0092 
0093   //  if(trhv->at(triggerNames_.triggerIndex("HLT_Mu9")).accept()) NHltMu9++;
0094   //  if(trhv->at(triggerNames_.triggerIndex("HLT_DoubleMu3")).accept()) NHltDiMu3++;
0095 
0096   this->handleObjects(iEvent, iSetup);
0097 
0098   if (useElectronSelection) {
0099     int nElecPassingCut = 0;
0100     for (unsigned int i = 0; i < theElectronCollection->size(); i++) {
0101       GsfElectron electron = (*theElectronCollection)[i];
0102       //      if (electron.ecalDrivenSeed()) {
0103       float elpt = electron.pt();
0104       //        if (electron.sigmaIetaIeta() < 0.002) continue;
0105       if (elpt > ptElecMin) {
0106         NtotalElectrons++;
0107         nElecPassingCut++;
0108         if (electron.mva_e_pi() > -0.1)
0109           NmvaElectrons++;
0110       }
0111       LogDebug("LeptonRecoSkim") << "elpt = " << elpt << endl;
0112       // } // closes if (electron.ecalDrivenSeed()) {
0113     }
0114     if (nElecPassingCut >= nSelectedElectrons)
0115       ElectronCutPassed = true;
0116   } else
0117     ElectronCutPassed = true;
0118 
0119   if (usePfElectronSelection) {
0120     int nPfElecPassingCut = 0;
0121     for (unsigned int i = 0; i < thePfCandidateCollection->size(); i++) {
0122       const reco::PFCandidate& thePfCandidate = (*thePfCandidateCollection)[i];
0123       if (thePfCandidate.particleId() != reco::PFCandidate::e)
0124         continue;
0125       if (thePfCandidate.gsfTrackRef().isNull())
0126         continue;
0127       float pfelpt = thePfCandidate.pt();
0128       //        if (electron.sigmaIetaIeta() < 0.002) continue;
0129       if (pfelpt > ptPfElecMin)
0130         nPfElecPassingCut++;
0131       LogDebug("LeptonRecoSkim") << "pfelpt = " << pfelpt << endl;
0132     }
0133     if (nPfElecPassingCut >= nSelectedPfElectrons)
0134       PfElectronCutPassed = true;
0135   } else
0136     PfElectronCutPassed = true;
0137 
0138   if (useMuonSelection) {
0139     int nMuonPassingCut = 0;
0140     for (unsigned int i = 0; i < theMuonCollection->size(); i++) {
0141       Muon muon = (*theMuonCollection)[i];
0142       if (!(muon.isGlobalMuon() || muon.isTrackerMuon()))
0143         continue;
0144       const TrackRef siTrack = muon.innerTrack();
0145       const TrackRef globalTrack = muon.globalTrack();
0146       float muonpt;
0147       float ptMuonMin;
0148       //      if (siTrack.isNonnull() && globalTrack.isNonnull()) {
0149       if (muon.isGlobalMuon() && muon.isTrackerMuon()) {
0150         muonpt = max(siTrack->pt(), globalTrack->pt());
0151         ptMuonMin = ptGlobalMuonMin;
0152       } else if (muon.isGlobalMuon() && !(muon.isTrackerMuon())) {
0153         muonpt = globalTrack->pt();
0154         ptMuonMin = ptGlobalMuonMin;
0155       } else if (muon.isTrackerMuon() && !(muon.isGlobalMuon())) {
0156         muonpt = siTrack->pt();
0157         ptMuonMin = ptTrackerMuonMin;
0158       } else {
0159         muonpt = 0;  // if we get here this is a STA only muon
0160         ptMuonMin = 999999.;
0161       }
0162       if (muonpt > ptMuonMin)
0163         nMuonPassingCut++;
0164       LogDebug("RecoSelectorCuts") << "muonpt = " << muonpt << endl;
0165     }
0166     if (nMuonPassingCut >= nSelectedMuons)
0167       MuonCutPassed = true;
0168   } else
0169     MuonCutPassed = true;
0170 
0171   if (useHtSelection) {
0172     double Ht = 0;
0173     for (unsigned int i = 0; i < theCaloJetCollection->size(); i++) {
0174       //      if((*theCaloJetCollection)[i].eta()<2.6 && (*theCaloJetCollection)[i].emEnergyFraction() <= 0.01) continue;
0175       if ((*theCaloJetCollection)[i].pt() > htJetThreshold)
0176         Ht += (*theCaloJetCollection)[i].pt();
0177     }
0178     if (Ht > htMin)
0179       HtCutPassed = true;
0180   } else
0181     HtCutPassed = true;
0182 
0183   if (usePFHtSelection) {
0184     double PFHt = 0;
0185     for (unsigned int i = 0; i < thePFJetCollection->size(); i++) {
0186       if ((*thePFJetCollection)[i].pt() > pfHtJetThreshold)
0187         PFHt += (*thePFJetCollection)[i].pt();
0188     }
0189     if (PFHt > pfHtMin)
0190       PFHtCutPassed = true;
0191   } else
0192     PFHtCutPassed = true;
0193 
0194   if (PfElectronCutPassed && ElectronCutPassed && MuonCutPassed && HtCutPassed && PFHtCutPassed)
0195     accept = true;
0196 
0197   if (accept)
0198     NeventsFiltered++;
0199 
0200   firstEvent = false;
0201   return accept;
0202 }
0203 
0204 // ------------ method called once each job just before starting event loop  ------------
0205 void LeptonRecoSkim::beginJob() { firstEvent = true; }
0206 
0207 // ------------ method called once each job just after ending the event loop  ------------
0208 void LeptonRecoSkim::endJob() {
0209   edm::LogPrint("LeptonRecoSkim") << "Filter Name            = " << filterName << endl;
0210   edm::LogPrint("LeptonRecoSkim") << "Total number of events = " << NeventsTotal << endl;
0211   edm::LogPrint("LeptonRecoSkim") << "Total HLT_Mu9          = " << NHltMu9 << endl;
0212   edm::LogPrint("LeptonRecoSkim") << "Total HLT_DoubleMu3    = " << NHltDiMu3 << endl;
0213   edm::LogPrint("LeptonRecoSkim") << "Filtered events        = " << NeventsFiltered << endl;
0214   edm::LogPrint("LeptonRecoSkim") << "Filter Efficiency      = " << (float)NeventsFiltered / (float)NeventsTotal
0215                                   << endl;
0216   edm::LogPrint("LeptonRecoSkim") << endl;
0217   edm::LogPrint("LeptonRecoSkim") << "N total electrons      = " << NtotalElectrons << endl;
0218   edm::LogPrint("LeptonRecoSkim") << "N mva>-0.1 electrons   = " << NmvaElectrons << endl;
0219   edm::LogPrint("LeptonRecoSkim") << endl;
0220   edm::LogPrint("LeptonRecoSkim") << endl;
0221 }
0222 
0223 void LeptonRecoSkim::handleObjects(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0224   //Get the electrons
0225   theElectronCollection = &iEvent.get(gsfElectronCollectionToken_);
0226 
0227   //Get the pf electrons
0228   thePfCandidateCollection = &iEvent.get(pfCandidateCollectionToken_);
0229 
0230   //Get the Muons
0231   theMuonCollection = &iEvent.get(muonCollectionToken_);
0232 
0233   //Get the CaloJets
0234   theCaloJetCollection = &iEvent.get(caloJetCollectionToken_);
0235 
0236   //Get the PfJets
0237   thePFJetCollection = &iEvent.get(pfJetCollectionToken_);
0238 
0239   // Get the ECAL rechhits to clean the spikes
0240   // Get EB RecHits
0241   theEcalBarrelCollection = &iEvent.get(ebRecHitCollectionToken_);
0242   // Get EE RecHits
0243   theEcalEndcapCollection = &iEvent.get(eeRecHitCollectionToken_);
0244 
0245   // Get topology for spike cleaning
0246   theCaloGeometry = &iSetup.getData(m_CaloGeoToken);
0247   theCaloTopology = &iSetup.getData(m_CaloTopoToken);
0248 }
0249 
0250 //define this as a plug-in
0251 DEFINE_FWK_MODULE(LeptonRecoSkim);